6d4a328b152dd935b8560579c0541ceecee49ccd
[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 "gromacs/timing/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
54 /* Bounding box calculations are (currently) always in single precision.
55  * This uses less (cache-)memory and SIMD is faster, at least on x86.
56  */
57 #define GMX_SIMD4_SINGLE
58 /* Include the 4-wide SIMD macro file */
59 #include "gmx_simd4_macros.h"
60 /* Check if we have 4-wide SIMD macro support */
61 #ifdef GMX_HAVE_SIMD4_MACROS
62 #define NBNXN_SEARCH_BB_SIMD4
63 #endif
64
65
66 #ifdef __cplusplus
67 extern "C" {
68 #endif
69
70
71 #ifdef GMX_NBNXN_SIMD
72 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
73 #define NBNXN_MEM_ALIGN  (GMX_SIMD_WIDTH_HERE*sizeof(real))
74 #else
75 /* No alignment required, but set it so we can call the same routines */
76 #define NBNXN_MEM_ALIGN  32
77 #endif
78
79
80 #ifdef NBNXN_SEARCH_BB_SIMD4
81 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
82 #define NBNXN_SEARCH_BB_MEM_ALIGN  (GMX_SIMD4_WIDTH*sizeof(float))
83 #else
84 /* No alignment required, but set it so we can call the same routines */
85 #define NBNXN_SEARCH_BB_MEM_ALIGN  32
86 #endif
87
88
89 /* Pair search box lower and upper corner in x,y,z.
90  * Store this in 4 iso 3 reals, which is useful with 4-wide SIMD.
91  * To avoid complicating the code we also use 4 without 4-wide SIMD.
92  */
93 #define NNBSBB_C         4
94 /* Pair search box lower and upper bound in z only. */
95 #define NNBSBB_D         2
96 /* Pair search box lower and upper corner x,y,z indices, entry 3 is unused */
97 #define BB_X  0
98 #define BB_Y  1
99 #define BB_Z  2
100
101 /* Bounding box for a nbnxn atom cluster */
102 typedef struct {
103     float lower[NNBSBB_C];
104     float upper[NNBSBB_C];
105 } nbnxn_bb_t;
106
107
108 /* A pair-search grid struct for one domain decomposition zone */
109 typedef struct {
110     rvec     c0;               /* The lower corner of the (local) grid        */
111     rvec     c1;               /* The upper corner of the (local) grid        */
112     real     atom_density;     /* The atom number density for the local grid  */
113
114     gmx_bool bSimple;          /* Is this grid simple or super/sub            */
115     int      na_c;             /* Number of atoms per cluster                 */
116     int      na_cj;            /* Number of atoms for list j-clusters         */
117     int      na_sc;            /* Number of atoms per super-cluster           */
118     int      na_c_2log;        /* 2log of na_c                                */
119
120     int      ncx;              /* Number of (super-)cells along x             */
121     int      ncy;              /* Number of (super-)cells along y             */
122     int      nc;               /* Total number of (super-)cells               */
123
124     real     sx;               /* x-size of a (super-)cell                    */
125     real     sy;               /* y-size of a (super-)cell                    */
126     real     inv_sx;           /* 1/sx                                        */
127     real     inv_sy;           /* 1/sy                                        */
128
129     int      cell0;            /* Index in nbs->cell corresponding to cell 0  */
130
131     int     *cxy_na;           /* The number of atoms for each column in x,y  */
132     int     *cxy_ind;          /* Grid (super)cell index, offset from cell0   */
133     int      cxy_nalloc;       /* Allocation size for cxy_na and cxy_ind      */
134
135     int        *nsubc;         /* The number of sub cells for each super cell */
136     float      *bbcz;          /* Bounding boxes in z for the super cells     */
137     nbnxn_bb_t *bb;            /* 3D bounding boxes for the sub cells         */
138     nbnxn_bb_t *bbj;           /* 3D j-bounding boxes for the case where      *
139                                 * the i- and j-cluster sizes are different    */
140     float      *pbb;           /* 3D b. boxes in xxxx format per super cell   */
141     int        *flags;         /* Flag for the super cells                    */
142     int         nc_nalloc;     /* Allocation size for the pointers above      */
143
144     float      *bbcz_simple;   /* bbcz for simple grid converted from super   */
145     nbnxn_bb_t *bb_simple;     /* bb for simple grid converted from super     */
146     int        *flags_simple;  /* flags for simple grid converted from super  */
147     int         nc_nalloc_simple; /* Allocation size for the pointers above   */
148
149     int      nsubc_tot;        /* Total number of subcell, used for printing  */
150 } nbnxn_grid_t;
151
152 #ifdef GMX_NBNXN_SIMD
153
154 typedef struct nbnxn_x_ci_simd_4xn {
155     /* The i-cluster coordinates for simple search */
156     gmx_mm_pr ix_S0, iy_S0, iz_S0;
157     gmx_mm_pr ix_S1, iy_S1, iz_S1;
158     gmx_mm_pr ix_S2, iy_S2, iz_S2;
159     gmx_mm_pr ix_S3, iy_S3, iz_S3;
160 } nbnxn_x_ci_simd_4xn_t;
161
162 typedef struct nbnxn_x_ci_simd_2xnn {
163     /* The i-cluster coordinates for simple search */
164     gmx_mm_pr ix_S0, iy_S0, iz_S0;
165     gmx_mm_pr ix_S2, iy_S2, iz_S2;
166 } nbnxn_x_ci_simd_2xnn_t;
167
168 #endif
169
170 /* Working data for the actual i-supercell during pair search */
171 typedef struct nbnxn_list_work {
172     gmx_cache_protect_t     cp0;    /* Protect cache between threads               */
173
174     nbnxn_bb_t             *bb_ci;  /* The bounding boxes, pbc shifted, for each cluster */
175     float                  *pbb_ci; /* As bb_ci, but in xxxx packed format               */
176     real                   *x_ci;   /* The coordinates, pbc shifted, for each atom       */
177 #ifdef GMX_NBNXN_SIMD
178     nbnxn_x_ci_simd_4xn_t  *x_ci_simd_4xn;
179     nbnxn_x_ci_simd_2xnn_t *x_ci_simd_2xnn;
180 #endif
181     int                     cj_ind;          /* The current cj_ind index for the current list     */
182     int                     cj4_init;        /* The first unitialized cj4 block                   */
183
184     float                  *d2;              /* Bounding box distance work array                  */
185
186     nbnxn_cj_t             *cj;              /* The j-cell list                                   */
187     int                     cj_nalloc;       /* Allocation size of cj                             */
188
189     int                     ncj_noq;         /* Nr. of cluster pairs without Coul for flop count  */
190     int                     ncj_hlj;         /* Nr. of cluster pairs with 1/2 LJ for flop count   */
191
192     int                    *sort;            /* Sort index                    */
193     int                     sort_nalloc;     /* Allocation size of sort       */
194
195     nbnxn_sci_t            *sci_sort;        /* Second sci array, for sorting */
196     int                     sci_sort_nalloc; /* Allocation size of sci_sort   */
197
198     gmx_cache_protect_t     cp1;             /* Protect cache between threads               */
199 } nbnxn_list_work_t;
200
201 /* Function type for setting the i-atom coordinate working data */
202 typedef void
203     gmx_icell_set_x_t (int ci,
204                        real shx, real shy, real shz,
205                        int na_c,
206                        int stride, const real *x,
207                        nbnxn_list_work_t *work);
208
209 static gmx_icell_set_x_t icell_set_x_simple;
210 #ifdef GMX_NBNXN_SIMD
211 static gmx_icell_set_x_t icell_set_x_simple_simd_4xn;
212 static gmx_icell_set_x_t icell_set_x_simple_simd_2xnn;
213 #endif
214 static gmx_icell_set_x_t icell_set_x_supersub;
215 #ifdef NBNXN_SEARCH_SSE
216 static gmx_icell_set_x_t icell_set_x_supersub_sse8;
217 #endif
218
219 /* Local cycle count struct for profiling */
220 typedef struct {
221     int          count;
222     gmx_cycles_t c;
223     gmx_cycles_t start;
224 } nbnxn_cycle_t;
225
226 /* Local cycle count enum for profiling */
227 enum {
228     enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
229 };
230
231 /* Thread-local work struct, contains part of nbnxn_grid_t */
232 typedef struct {
233     gmx_cache_protect_t  cp0;
234
235     int                 *cxy_na;
236     int                  cxy_na_nalloc;
237
238     int                 *sort_work;
239     int                  sort_work_nalloc;
240
241     nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
242
243     int                  ndistc;       /* Number of distance checks for flop counting */
244
245     nbnxn_cycle_t        cc[enbsCCnr];
246
247     gmx_cache_protect_t  cp1;
248 } nbnxn_search_work_t;
249
250 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
251 typedef struct nbnxn_search {
252     int                 ePBC;            /* PBC type enum                              */
253     matrix              box;             /* The periodic unit-cell                     */
254
255     gmx_bool            DomDec;          /* Are we doing domain decomposition?         */
256     ivec                dd_dim;          /* Are we doing DD in x,y,z?                  */
257     gmx_domdec_zones_t *zones;           /* The domain decomposition zones        */
258
259     int                 ngrid;           /* The number of grids, equal to #DD-zones    */
260     nbnxn_grid_t       *grid;            /* Array of grids, size ngrid                 */
261     int                *cell;            /* Actual allocated cell array for all grids  */
262     int                 cell_nalloc;     /* Allocation size of cell                    */
263     int                *a;               /* Atom index for grid, the inverse of cell   */
264     int                 a_nalloc;        /* Allocation size of a                       */
265
266     int                 natoms_local;    /* The local atoms run from 0 to natoms_local */
267     int                 natoms_nonlocal; /* The non-local atoms run from natoms_local
268                                           * to natoms_nonlocal */
269
270     gmx_bool             print_cycles;
271     int                  search_count;
272     nbnxn_cycle_t        cc[enbsCCnr];
273
274     gmx_icell_set_x_t   *icell_set_x; /* Function for setting i-coords    */
275
276     int                  nthread_max; /* Maximum number of threads for pair-search  */
277     nbnxn_search_work_t *work;        /* Work array, size nthread_max          */
278 } nbnxn_search_t_t;
279
280
281 static void nbs_cycle_start(nbnxn_cycle_t *cc)
282 {
283     cc->start = gmx_cycles_read();
284 }
285
286 static void nbs_cycle_stop(nbnxn_cycle_t *cc)
287 {
288     cc->c += gmx_cycles_read() - cc->start;
289     cc->count++;
290 }
291
292
293 #ifdef __cplusplus
294 }
295 #endif
296
297 #endif