Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_internal.h
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustr
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2012, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14  *
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35
36 #ifndef _nbnxn_internal_h
37 #define _nbnxn_internal_h
38
39 #include "typedefs.h"
40 #include "domdec.h"
41 #include "gmx_cyclecounter.h"
42
43 #ifdef GMX_NBNXN_SIMD
44 /* The include below sets the SIMD instruction type (precision+width)
45  * for all nbnxn SIMD search and non-bonded kernel code.
46  */
47 #ifdef GMX_NBNXN_HALF_WIDTH_SIMD
48 #define GMX_USE_HALF_WIDTH_SIMD_HERE
49 #endif
50 #include "gmx_simd_macros.h"
51 #endif
52
53 #ifdef __cplusplus
54 extern "C" {
55 #endif
56
57
58 #ifdef GMX_X86_SSE2
59 /* Use 4-way SIMD for, always, single precision bounding box calculations */
60 #define NBNXN_SEARCH_BB_SSE
61 #endif
62
63
64 #ifdef GMX_NBNXN_SIMD
65 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
66 #define NBNXN_MEM_ALIGN  (GMX_SIMD_WIDTH_HERE*sizeof(real))
67 #else
68 /* No alignment required, but set it so we can call the same routines */
69 #define NBNXN_MEM_ALIGN  32
70 #endif
71
72
73 /* Pair search box lower and upper corner in x,y,z.
74  * Store this in 4 iso 3 reals, which is useful with SSE.
75  * To avoid complicating the code we also use 4 without SSE.
76  */
77 #define NNBSBB_C         4
78 /* Pair search box lower and upper bound in z only. */
79 #define NNBSBB_D         2
80 /* Pair search box lower and upper corner x,y,z indices, entry 3 is unused */
81 #define BB_X  0
82 #define BB_Y  1
83 #define BB_Z  2
84
85 /* Bounding box for a nbnxn atom cluster */
86 typedef struct {
87     float lower[NNBSBB_C];
88     float upper[NNBSBB_C];
89 } nbnxn_bb_t;
90
91
92 /* A pair-search grid struct for one domain decomposition zone */
93 typedef struct {
94     rvec     c0;               /* The lower corner of the (local) grid        */
95     rvec     c1;               /* The upper corner of the (local) grid        */
96     real     atom_density;     /* The atom number density for the local grid  */
97
98     gmx_bool bSimple;          /* Is this grid simple or super/sub            */
99     int      na_c;             /* Number of atoms per cluster                 */
100     int      na_cj;            /* Number of atoms for list j-clusters         */
101     int      na_sc;            /* Number of atoms per super-cluster           */
102     int      na_c_2log;        /* 2log of na_c                                */
103
104     int      ncx;              /* Number of (super-)cells along x             */
105     int      ncy;              /* Number of (super-)cells along y             */
106     int      nc;               /* Total number of (super-)cells               */
107
108     real     sx;               /* x-size of a (super-)cell                    */
109     real     sy;               /* y-size of a (super-)cell                    */
110     real     inv_sx;           /* 1/sx                                        */
111     real     inv_sy;           /* 1/sy                                        */
112
113     int      cell0;            /* Index in nbs->cell corresponding to cell 0  */
114
115     int     *cxy_na;           /* The number of atoms for each column in x,y  */
116     int     *cxy_ind;          /* Grid (super)cell index, offset from cell0   */
117     int      cxy_nalloc;       /* Allocation size for cxy_na and cxy_ind      */
118
119     int        *nsubc;         /* The number of sub cells for each super cell */
120     float      *bbcz;          /* Bounding boxes in z for the super cells     */
121     nbnxn_bb_t *bb;            /* 3D bounding boxes for the sub cells         */
122     nbnxn_bb_t *bbj;           /* 3D j-b.boxes for SSE-double or AVX-single   */
123     float      *pbb;           /* 3D b. boxes in xxxx format per super cell   */
124     int        *flags;         /* Flag for the super cells                    */
125     int         nc_nalloc;     /* Allocation size for the pointers above      */
126
127     float      *bbcz_simple;   /* bbcz for simple grid converted from super   */
128     nbnxn_bb_t *bb_simple;     /* bb for simple grid converted from super     */
129     int        *flags_simple;  /* flags for simple grid converted from super  */
130     int         nc_nalloc_simple; /* Allocation size for the pointers above   */
131
132     int      nsubc_tot;        /* Total number of subcell, used for printing  */
133 } nbnxn_grid_t;
134
135 #ifdef GMX_NBNXN_SIMD
136
137 typedef struct nbnxn_x_ci_simd_4xn {
138     /* The i-cluster coordinates for simple search */
139     gmx_mm_pr ix_SSE0, iy_SSE0, iz_SSE0;
140     gmx_mm_pr ix_SSE1, iy_SSE1, iz_SSE1;
141     gmx_mm_pr ix_SSE2, iy_SSE2, iz_SSE2;
142     gmx_mm_pr ix_SSE3, iy_SSE3, iz_SSE3;
143 } nbnxn_x_ci_simd_4xn_t;
144
145 typedef struct nbnxn_x_ci_simd_2xnn {
146     /* The i-cluster coordinates for simple search */
147     gmx_mm_pr ix_SSE0, iy_SSE0, iz_SSE0;
148     gmx_mm_pr ix_SSE2, iy_SSE2, iz_SSE2;
149 } nbnxn_x_ci_simd_2xnn_t;
150
151 #endif
152
153 /* Working data for the actual i-supercell during pair search */
154 typedef struct nbnxn_list_work {
155     gmx_cache_protect_t     cp0;    /* Protect cache between threads               */
156
157     nbnxn_bb_t             *bb_ci;  /* The bounding boxes, pbc shifted, for each cluster */
158     float                  *pbb_ci; /* As bb_ci, but in xxxx packed format               */
159     real                   *x_ci;   /* The coordinates, pbc shifted, for each atom       */
160 #ifdef GMX_NBNXN_SIMD
161     nbnxn_x_ci_simd_4xn_t  *x_ci_simd_4xn;
162     nbnxn_x_ci_simd_2xnn_t *x_ci_simd_2xnn;
163 #endif
164     int                     cj_ind;          /* The current cj_ind index for the current list     */
165     int                     cj4_init;        /* The first unitialized cj4 block                   */
166
167     float                  *d2;              /* Bounding box distance work array                  */
168
169     nbnxn_cj_t             *cj;              /* The j-cell list                                   */
170     int                     cj_nalloc;       /* Allocation size of cj                             */
171
172     int                     ncj_noq;         /* Nr. of cluster pairs without Coul for flop count  */
173     int                     ncj_hlj;         /* Nr. of cluster pairs with 1/2 LJ for flop count   */
174
175     int                    *sort;            /* Sort index                    */
176     int                     sort_nalloc;     /* Allocation size of sort       */
177
178     nbnxn_sci_t            *sci_sort;        /* Second sci array, for sorting */
179     int                     sci_sort_nalloc; /* Allocation size of sci_sort   */
180
181     gmx_cache_protect_t     cp1;             /* Protect cache between threads               */
182 } nbnxn_list_work_t;
183
184 /* Function type for setting the i-atom coordinate working data */
185 typedef void
186     gmx_icell_set_x_t (int ci,
187                        real shx, real shy, real shz,
188                        int na_c,
189                        int stride, const real *x,
190                        nbnxn_list_work_t *work);
191
192 static gmx_icell_set_x_t icell_set_x_simple;
193 #ifdef GMX_NBNXN_SIMD
194 static gmx_icell_set_x_t icell_set_x_simple_simd_4xn;
195 static gmx_icell_set_x_t icell_set_x_simple_simd_2xnn;
196 #endif
197 static gmx_icell_set_x_t icell_set_x_supersub;
198 #ifdef NBNXN_SEARCH_SSE
199 static gmx_icell_set_x_t icell_set_x_supersub_sse8;
200 #endif
201
202 /* Local cycle count struct for profiling */
203 typedef struct {
204     int          count;
205     gmx_cycles_t c;
206     gmx_cycles_t start;
207 } nbnxn_cycle_t;
208
209 /* Local cycle count enum for profiling */
210 enum {
211     enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
212 };
213
214 /* Thread-local work struct, contains part of nbnxn_grid_t */
215 typedef struct {
216     gmx_cache_protect_t  cp0;
217
218     int                 *cxy_na;
219     int                  cxy_na_nalloc;
220
221     int                 *sort_work;
222     int                  sort_work_nalloc;
223
224     nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
225
226     int                  ndistc;       /* Number of distance checks for flop counting */
227
228     nbnxn_cycle_t        cc[enbsCCnr];
229
230     gmx_cache_protect_t  cp1;
231 } nbnxn_search_work_t;
232
233 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
234 typedef struct nbnxn_search {
235     int                 ePBC;            /* PBC type enum                              */
236     matrix              box;             /* The periodic unit-cell                     */
237
238     gmx_bool            DomDec;          /* Are we doing domain decomposition?         */
239     ivec                dd_dim;          /* Are we doing DD in x,y,z?                  */
240     gmx_domdec_zones_t *zones;           /* The domain decomposition zones        */
241
242     int                 ngrid;           /* The number of grids, equal to #DD-zones    */
243     nbnxn_grid_t       *grid;            /* Array of grids, size ngrid                 */
244     int                *cell;            /* Actual allocated cell array for all grids  */
245     int                 cell_nalloc;     /* Allocation size of cell                    */
246     int                *a;               /* Atom index for grid, the inverse of cell   */
247     int                 a_nalloc;        /* Allocation size of a                       */
248
249     int                 natoms_local;    /* The local atoms run from 0 to natoms_local */
250     int                 natoms_nonlocal; /* The non-local atoms run from natoms_local
251                                           * to natoms_nonlocal */
252
253     gmx_bool             print_cycles;
254     int                  search_count;
255     nbnxn_cycle_t        cc[enbsCCnr];
256
257     gmx_icell_set_x_t   *icell_set_x; /* Function for setting i-coords    */
258
259     int                  nthread_max; /* Maximum number of threads for pair-search  */
260     nbnxn_search_work_t *work;        /* Work array, size nthread_max          */
261 } nbnxn_search_t_t;
262
263
264 static void nbs_cycle_start(nbnxn_cycle_t *cc)
265 {
266     cc->start = gmx_cycles_read();
267 }
268
269 static void nbs_cycle_stop(nbnxn_cycle_t *cc)
270 {
271     cc->c += gmx_cycles_read() - cc->start;
272     cc->count++;
273 }
274
275
276 #ifdef __cplusplus
277 }
278 #endif
279
280 #endif