5d5844cb226bd69ed37b454fe04c4387880502ea
[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/legacyheaders/typedefs.h"
40 #include "nbnxn_simd.h"
41 #include "gromacs/legacyheaders/domdec.h"
42 #include "gromacs/timing/cyclecounter.h"
43 #include "gromacs/mdlib/nbnxn_pairlist.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     real          atom_density;     /* The atom number density for the local grid  */
102
103     gmx_bool      bSimple;          /* Is this grid simple or super/sub            */
104     int           na_c;             /* Number of atoms per cluster                 */
105     int           na_cj;            /* Number of atoms for list j-clusters         */
106     int           na_sc;            /* Number of atoms per super-cluster           */
107     int           na_c_2log;        /* 2log of na_c                                */
108
109     int           ncx;              /* Number of (super-)cells along x             */
110     int           ncy;              /* Number of (super-)cells along y             */
111     int           nc;               /* Total number of (super-)cells               */
112
113     real          sx;               /* x-size of a (super-)cell                    */
114     real          sy;               /* y-size of a (super-)cell                    */
115     real          inv_sx;           /* 1/sx                                        */
116     real          inv_sy;           /* 1/sy                                        */
117
118     int           cell0;            /* Index in nbs->cell corresponding to cell 0  */
119
120     int          *cxy_na;           /* The number of atoms for each column in x,y  */
121     int          *cxy_ind;          /* Grid (super)cell index, offset from cell0   */
122     int           cxy_nalloc;       /* Allocation size for cxy_na and cxy_ind      */
123
124     int          *nsubc;            /* The number of sub cells for each super cell */
125     float        *bbcz;             /* Bounding boxes in z for the super cells     */
126     nbnxn_bb_t   *bb;               /* 3D bounding boxes for the sub cells         */
127     nbnxn_bb_t   *bbj;              /* 3D j-bounding boxes for the case where      *
128                                      * the i- and j-cluster sizes are different    */
129     float        *pbb;              /* 3D b. boxes in xxxx format per super cell   */
130     int          *flags;            /* Flag for the super cells                    */
131     unsigned int *fep;              /* FEP signal bits for sub cells               */
132     int           nc_nalloc;        /* Allocation size for the pointers above      */
133
134     float        *bbcz_simple;      /* bbcz for simple grid converted from super   */
135     nbnxn_bb_t   *bb_simple;        /* bb for simple grid converted from super     */
136     int          *flags_simple;     /* flags for simple grid converted from super  */
137     int           nc_nalloc_simple; /* Allocation size for the pointers above   */
138
139     int           nsubc_tot;        /* Total number of subcell, used for printing  */
140 } nbnxn_grid_t;
141
142 #ifdef GMX_NBNXN_SIMD
143
144 typedef struct nbnxn_x_ci_simd_4xn {
145     /* The i-cluster coordinates for simple search */
146     gmx_simd_real_t ix_S0, iy_S0, iz_S0;
147     gmx_simd_real_t ix_S1, iy_S1, iz_S1;
148     gmx_simd_real_t ix_S2, iy_S2, iz_S2;
149     gmx_simd_real_t ix_S3, iy_S3, iz_S3;
150 } nbnxn_x_ci_simd_4xn_t;
151
152 typedef struct nbnxn_x_ci_simd_2xnn {
153     /* The i-cluster coordinates for simple search */
154     gmx_simd_real_t ix_S0, iy_S0, iz_S0;
155     gmx_simd_real_t ix_S2, iy_S2, iz_S2;
156 } nbnxn_x_ci_simd_2xnn_t;
157
158 #endif
159
160 /* Working data for the actual i-supercell during pair search */
161 typedef struct nbnxn_list_work {
162     gmx_cache_protect_t     cp0;    /* Protect cache between threads               */
163
164     nbnxn_bb_t             *bb_ci;  /* The bounding boxes, pbc shifted, for each cluster */
165     float                  *pbb_ci; /* As bb_ci, but in xxxx packed format               */
166     real                   *x_ci;   /* The coordinates, pbc shifted, for each atom       */
167 #ifdef GMX_NBNXN_SIMD
168     nbnxn_x_ci_simd_4xn_t  *x_ci_simd_4xn;
169     nbnxn_x_ci_simd_2xnn_t *x_ci_simd_2xnn;
170 #endif
171     int                     cj_ind;          /* The current cj_ind index for the current list     */
172     int                     cj4_init;        /* The first unitialized cj4 block                   */
173
174     float                  *d2;              /* Bounding box distance work array                  */
175
176     nbnxn_cj_t             *cj;              /* The j-cell list                                   */
177     int                     cj_nalloc;       /* Allocation size of cj                             */
178
179     int                     ncj_noq;         /* Nr. of cluster pairs without Coul for flop count  */
180     int                     ncj_hlj;         /* Nr. of cluster pairs with 1/2 LJ for flop count   */
181
182     int                    *sort;            /* Sort index                    */
183     int                     sort_nalloc;     /* Allocation size of sort       */
184
185     nbnxn_sci_t            *sci_sort;        /* Second sci array, for sorting */
186     int                     sci_sort_nalloc; /* Allocation size of sci_sort   */
187
188     gmx_cache_protect_t     cp1;             /* Protect cache between threads               */
189 } nbnxn_list_work_t;
190
191 /* Function type for setting the i-atom coordinate working data */
192 typedef void
193     gmx_icell_set_x_t (int ci,
194                        real shx, real shy, real shz,
195                        int na_c,
196                        int stride, const real *x,
197                        nbnxn_list_work_t *work);
198
199 static gmx_icell_set_x_t icell_set_x_simple;
200 #ifdef GMX_NBNXN_SIMD
201 static gmx_icell_set_x_t icell_set_x_simple_simd_4xn;
202 static gmx_icell_set_x_t icell_set_x_simple_simd_2xnn;
203 #endif
204 static gmx_icell_set_x_t icell_set_x_supersub;
205 #ifdef NBNXN_SEARCH_SSE
206 static gmx_icell_set_x_t icell_set_x_supersub_sse8;
207 #endif
208
209 /* Local cycle count struct for profiling */
210 typedef struct {
211     int          count;
212     gmx_cycles_t c;
213     gmx_cycles_t start;
214 } nbnxn_cycle_t;
215
216 /* Local cycle count enum for profiling */
217 enum {
218     enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
219 };
220
221 /* Thread-local work struct, contains part of nbnxn_grid_t */
222 typedef struct {
223     gmx_cache_protect_t  cp0;
224
225     int                 *cxy_na;
226     int                  cxy_na_nalloc;
227
228     int                 *sort_work;
229     int                  sort_work_nalloc;
230
231     nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
232
233     int                  ndistc;       /* Number of distance checks for flop counting */
234
235     t_nblist            *nbl_fep;      /* Temporary FEP list for load balancing */
236
237     nbnxn_cycle_t        cc[enbsCCnr];
238
239     gmx_cache_protect_t  cp1;
240 } nbnxn_search_work_t;
241
242 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
243 typedef struct nbnxn_search {
244     gmx_bool            bFEP;            /* Do we have perturbed atoms? */
245     int                 ePBC;            /* PBC type enum                              */
246     matrix              box;             /* The periodic unit-cell                     */
247
248     gmx_bool            DomDec;          /* Are we doing domain decomposition?         */
249     ivec                dd_dim;          /* Are we doing DD in x,y,z?                  */
250     gmx_domdec_zones_t *zones;           /* The domain decomposition zones        */
251
252     int                 ngrid;           /* The number of grids, equal to #DD-zones    */
253     nbnxn_grid_t       *grid;            /* Array of grids, size ngrid                 */
254     int                *cell;            /* Actual allocated cell array for all grids  */
255     int                 cell_nalloc;     /* Allocation size of cell                    */
256     int                *a;               /* Atom index for grid, the inverse of cell   */
257     int                 a_nalloc;        /* Allocation size of a                       */
258
259     int                 natoms_local;    /* The local atoms run from 0 to natoms_local */
260     int                 natoms_nonlocal; /* The non-local atoms run from natoms_local
261                                           * to natoms_nonlocal */
262
263     gmx_bool             print_cycles;
264     int                  search_count;
265     nbnxn_cycle_t        cc[enbsCCnr];
266
267     gmx_icell_set_x_t   *icell_set_x; /* Function for setting i-coords    */
268
269     int                  nthread_max; /* Maximum number of threads for pair-search  */
270     nbnxn_search_work_t *work;        /* Work array, size nthread_max          */
271 } nbnxn_search_t_t;
272
273
274 static void nbs_cycle_start(nbnxn_cycle_t *cc)
275 {
276     cc->start = gmx_cycles_read();
277 }
278
279 static void nbs_cycle_stop(nbnxn_cycle_t *cc)
280 {
281     cc->c += gmx_cycles_read() - cc->start;
282     cc->count++;
283 }
284
285
286 #ifdef __cplusplus
287 }
288 #endif
289
290 #endif