2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2012, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 #ifndef _nbnxn_internal_h
40 #define _nbnxn_internal_h
44 #include "gmx_cyclecounter.h"
47 /* The include below sets the SIMD instruction type (precision+width)
48 * for all nbnxn SIMD search and non-bonded kernel code.
50 #ifdef GMX_NBNXN_HALF_WIDTH_SIMD
51 #define GMX_USE_HALF_WIDTH_SIMD_HERE
53 #include "gmx_simd_macros.h"
57 /* Bounding box calculations are (currently) always in single precision.
58 * This uses less (cache-)memory and SIMD is faster, at least on x86.
60 #define GMX_SIMD4_SINGLE
61 /* Include the 4-wide SIMD macro file */
62 #include "gmx_simd4_macros.h"
63 /* Check if we have 4-wide SIMD macro support */
64 #ifdef GMX_HAVE_SIMD4_MACROS
65 #define NBNXN_SEARCH_BB_SIMD4
75 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
76 #define NBNXN_MEM_ALIGN (GMX_SIMD_WIDTH_HERE*sizeof(real))
78 /* No alignment required, but set it so we can call the same routines */
79 #define NBNXN_MEM_ALIGN 32
83 #ifdef NBNXN_SEARCH_BB_SIMD4
84 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
85 #define NBNXN_SEARCH_BB_MEM_ALIGN (GMX_SIMD4_WIDTH*sizeof(float))
87 /* No alignment required, but set it so we can call the same routines */
88 #define NBNXN_SEARCH_BB_MEM_ALIGN 32
92 /* Pair search box lower and upper corner in x,y,z.
93 * Store this in 4 iso 3 reals, which is useful with 4-wide SIMD.
94 * To avoid complicating the code we also use 4 without 4-wide SIMD.
97 /* Pair search box lower and upper bound in z only. */
99 /* Pair search box lower and upper corner x,y,z indices, entry 3 is unused */
104 /* Bounding box for a nbnxn atom cluster */
106 float lower[NNBSBB_C];
107 float upper[NNBSBB_C];
111 /* A pair-search grid struct for one domain decomposition zone */
113 rvec c0; /* The lower corner of the (local) grid */
114 rvec c1; /* The upper corner of the (local) grid */
115 real atom_density; /* The atom number density for the local grid */
117 gmx_bool bSimple; /* Is this grid simple or super/sub */
118 int na_c; /* Number of atoms per cluster */
119 int na_cj; /* Number of atoms for list j-clusters */
120 int na_sc; /* Number of atoms per super-cluster */
121 int na_c_2log; /* 2log of na_c */
123 int ncx; /* Number of (super-)cells along x */
124 int ncy; /* Number of (super-)cells along y */
125 int nc; /* Total number of (super-)cells */
127 real sx; /* x-size of a (super-)cell */
128 real sy; /* y-size of a (super-)cell */
129 real inv_sx; /* 1/sx */
130 real inv_sy; /* 1/sy */
132 int cell0; /* Index in nbs->cell corresponding to cell 0 */
134 int *cxy_na; /* The number of atoms for each column in x,y */
135 int *cxy_ind; /* Grid (super)cell index, offset from cell0 */
136 int cxy_nalloc; /* Allocation size for cxy_na and cxy_ind */
138 int *nsubc; /* The number of sub cells for each super cell */
139 float *bbcz; /* Bounding boxes in z for the super cells */
140 nbnxn_bb_t *bb; /* 3D bounding boxes for the sub cells */
141 nbnxn_bb_t *bbj; /* 3D j-bounding boxes for the case where *
142 * the i- and j-cluster sizes are different */
143 float *pbb; /* 3D b. boxes in xxxx format per super cell */
144 int *flags; /* Flag for the super cells */
145 int nc_nalloc; /* Allocation size for the pointers above */
147 float *bbcz_simple; /* bbcz for simple grid converted from super */
148 nbnxn_bb_t *bb_simple; /* bb for simple grid converted from super */
149 int *flags_simple; /* flags for simple grid converted from super */
150 int nc_nalloc_simple; /* Allocation size for the pointers above */
152 int nsubc_tot; /* Total number of subcell, used for printing */
155 #ifdef GMX_NBNXN_SIMD
157 typedef struct nbnxn_x_ci_simd_4xn {
158 /* The i-cluster coordinates for simple search */
159 gmx_mm_pr ix_S0, iy_S0, iz_S0;
160 gmx_mm_pr ix_S1, iy_S1, iz_S1;
161 gmx_mm_pr ix_S2, iy_S2, iz_S2;
162 gmx_mm_pr ix_S3, iy_S3, iz_S3;
163 } nbnxn_x_ci_simd_4xn_t;
165 typedef struct nbnxn_x_ci_simd_2xnn {
166 /* The i-cluster coordinates for simple search */
167 gmx_mm_pr ix_S0, iy_S0, iz_S0;
168 gmx_mm_pr ix_S2, iy_S2, iz_S2;
169 } nbnxn_x_ci_simd_2xnn_t;
173 /* Working data for the actual i-supercell during pair search */
174 typedef struct nbnxn_list_work {
175 gmx_cache_protect_t cp0; /* Protect cache between threads */
177 nbnxn_bb_t *bb_ci; /* The bounding boxes, pbc shifted, for each cluster */
178 float *pbb_ci; /* As bb_ci, but in xxxx packed format */
179 real *x_ci; /* The coordinates, pbc shifted, for each atom */
180 #ifdef GMX_NBNXN_SIMD
181 nbnxn_x_ci_simd_4xn_t *x_ci_simd_4xn;
182 nbnxn_x_ci_simd_2xnn_t *x_ci_simd_2xnn;
184 int cj_ind; /* The current cj_ind index for the current list */
185 int cj4_init; /* The first unitialized cj4 block */
187 float *d2; /* Bounding box distance work array */
189 nbnxn_cj_t *cj; /* The j-cell list */
190 int cj_nalloc; /* Allocation size of cj */
192 int ncj_noq; /* Nr. of cluster pairs without Coul for flop count */
193 int ncj_hlj; /* Nr. of cluster pairs with 1/2 LJ for flop count */
195 int *sort; /* Sort index */
196 int sort_nalloc; /* Allocation size of sort */
198 nbnxn_sci_t *sci_sort; /* Second sci array, for sorting */
199 int sci_sort_nalloc; /* Allocation size of sci_sort */
201 gmx_cache_protect_t cp1; /* Protect cache between threads */
204 /* Function type for setting the i-atom coordinate working data */
206 gmx_icell_set_x_t (int ci,
207 real shx, real shy, real shz,
209 int stride, const real *x,
210 nbnxn_list_work_t *work);
212 static gmx_icell_set_x_t icell_set_x_simple;
213 #ifdef GMX_NBNXN_SIMD
214 static gmx_icell_set_x_t icell_set_x_simple_simd_4xn;
215 static gmx_icell_set_x_t icell_set_x_simple_simd_2xnn;
217 static gmx_icell_set_x_t icell_set_x_supersub;
218 #ifdef NBNXN_SEARCH_SSE
219 static gmx_icell_set_x_t icell_set_x_supersub_sse8;
222 /* Local cycle count struct for profiling */
229 /* Local cycle count enum for profiling */
231 enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
234 /* Thread-local work struct, contains part of nbnxn_grid_t */
236 gmx_cache_protect_t cp0;
242 int sort_work_nalloc;
244 nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
246 int ndistc; /* Number of distance checks for flop counting */
248 nbnxn_cycle_t cc[enbsCCnr];
250 gmx_cache_protect_t cp1;
251 } nbnxn_search_work_t;
253 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
254 typedef struct nbnxn_search {
255 int ePBC; /* PBC type enum */
256 matrix box; /* The periodic unit-cell */
258 gmx_bool DomDec; /* Are we doing domain decomposition? */
259 ivec dd_dim; /* Are we doing DD in x,y,z? */
260 gmx_domdec_zones_t *zones; /* The domain decomposition zones */
262 int ngrid; /* The number of grids, equal to #DD-zones */
263 nbnxn_grid_t *grid; /* Array of grids, size ngrid */
264 int *cell; /* Actual allocated cell array for all grids */
265 int cell_nalloc; /* Allocation size of cell */
266 int *a; /* Atom index for grid, the inverse of cell */
267 int a_nalloc; /* Allocation size of a */
269 int natoms_local; /* The local atoms run from 0 to natoms_local */
270 int natoms_nonlocal; /* The non-local atoms run from natoms_local
271 * to natoms_nonlocal */
273 gmx_bool print_cycles;
275 nbnxn_cycle_t cc[enbsCCnr];
277 gmx_icell_set_x_t *icell_set_x; /* Function for setting i-coords */
279 int nthread_max; /* Maximum number of threads for pair-search */
280 nbnxn_search_work_t *work; /* Work array, size nthread_max */
284 static void nbs_cycle_start(nbnxn_cycle_t *cc)
286 cc->start = gmx_cycles_read();
289 static void nbs_cycle_stop(nbnxn_cycle_t *cc)
291 cc->c += gmx_cycles_read() - cc->start;