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"
52 /* Use 4-way SIMD for, always, single precision bounding box calculations */
53 #define NBNXN_SEARCH_BB_SSE
58 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
59 #define NBNXN_MEM_ALIGN (GMX_NBNXN_SIMD_BITWIDTH/8)
61 /* No alignment required, but set it so we can call the same routines */
62 #define NBNXN_MEM_ALIGN 32
66 /* A pair-search grid struct for one domain decomposition zone */
68 rvec c0; /* The lower corner of the (local) grid */
69 rvec c1; /* The upper corner of the (local) grid */
70 real atom_density; /* The atom number density for the local grid */
72 gmx_bool bSimple; /* Is this grid simple or super/sub */
73 int na_c; /* Number of atoms per cluster */
74 int na_cj; /* Number of atoms for list j-clusters */
75 int na_sc; /* Number of atoms per super-cluster */
76 int na_c_2log; /* 2log of na_c */
78 int ncx; /* Number of (super-)cells along x */
79 int ncy; /* Number of (super-)cells along y */
80 int nc; /* Total number of (super-)cells */
82 real sx; /* x-size of a (super-)cell */
83 real sy; /* y-size of a (super-)cell */
84 real inv_sx; /* 1/sx */
85 real inv_sy; /* 1/sy */
87 int cell0; /* Index in nbs->cell corresponding to cell 0 */
89 int *cxy_na; /* The number of atoms for each column in x,y */
90 int *cxy_ind; /* Grid (super)cell index, offset from cell0 */
91 int cxy_nalloc; /* Allocation size for cxy_na and cxy_ind */
93 int *nsubc; /* The number of sub cells for each super cell */
94 float *bbcz; /* Bounding boxes in z for the super cells */
95 float *bb; /* 3D bounding boxes for the sub cells */
96 float *bbj; /* 3D j-b.boxes for SSE-double or AVX-single */
97 int *flags; /* Flag for the super cells */
98 int nc_nalloc; /* Allocation size for the pointers above */
100 float *bbcz_simple; /* bbcz for simple grid converted from super */
101 float *bb_simple; /* bb for simple grid converted from super */
102 int *flags_simple; /* flags for simple grid converted from super */
103 int nc_nalloc_simple; /* Allocation size for the pointers above */
105 int nsubc_tot; /* Total number of subcell, used for printing */
108 #ifdef GMX_NBNXN_SIMD
109 #ifdef GMX_NBNXN_HALF_WIDTH_SIMD
110 #define GMX_USE_HALF_WIDTH_SIMD_HERE
112 #include "gmx_simd_macros.h"
114 typedef struct nbnxn_x_ci_simd_4xn {
115 /* The i-cluster coordinates for simple search */
116 gmx_mm_pr ix_SSE0, iy_SSE0, iz_SSE0;
117 gmx_mm_pr ix_SSE1, iy_SSE1, iz_SSE1;
118 gmx_mm_pr ix_SSE2, iy_SSE2, iz_SSE2;
119 gmx_mm_pr ix_SSE3, iy_SSE3, iz_SSE3;
120 } nbnxn_x_ci_simd_4xn_t;
122 typedef struct nbnxn_x_ci_simd_2xnn {
123 /* The i-cluster coordinates for simple search */
124 gmx_mm_pr ix_SSE0, iy_SSE0, iz_SSE0;
125 gmx_mm_pr ix_SSE2, iy_SSE2, iz_SSE2;
126 } nbnxn_x_ci_simd_2xnn_t;
130 /* Working data for the actual i-supercell during pair search */
131 typedef struct nbnxn_list_work {
132 gmx_cache_protect_t cp0; /* Protect cache between threads */
134 float *bb_ci; /* The bounding boxes, pbc shifted, for each cluster */
135 real *x_ci; /* The coordinates, pbc shifted, for each atom */
136 #ifdef GMX_NBNXN_SIMD
137 nbnxn_x_ci_simd_4xn_t *x_ci_simd_4xn;
138 nbnxn_x_ci_simd_2xnn_t *x_ci_simd_2xnn;
140 int cj_ind; /* The current cj_ind index for the current list */
141 int cj4_init; /* The first unitialized cj4 block */
143 float *d2; /* Bounding box distance work array */
145 nbnxn_cj_t *cj; /* The j-cell list */
146 int cj_nalloc; /* Allocation size of cj */
148 int ncj_noq; /* Nr. of cluster pairs without Coul for flop count */
149 int ncj_hlj; /* Nr. of cluster pairs with 1/2 LJ for flop count */
151 int *sort; /* Sort index */
152 int sort_nalloc; /* Allocation size of sort */
154 nbnxn_sci_t *sci_sort; /* Second sci array, for sorting */
155 int sci_sort_nalloc; /* Allocation size of sci_sort */
157 gmx_cache_protect_t cp1; /* Protect cache between threads */
160 /* Function type for setting the i-atom coordinate working data */
162 gmx_icell_set_x_t (int ci,
163 real shx, real shy, real shz,
165 int stride, const real *x,
166 nbnxn_list_work_t *work);
168 static gmx_icell_set_x_t icell_set_x_simple;
169 #ifdef GMX_NBNXN_SIMD
170 static gmx_icell_set_x_t icell_set_x_simple_simd_4xn;
171 static gmx_icell_set_x_t icell_set_x_simple_simd_2xnn;
173 static gmx_icell_set_x_t icell_set_x_supersub;
174 #ifdef NBNXN_SEARCH_SSE
175 static gmx_icell_set_x_t icell_set_x_supersub_sse8;
178 /* Local cycle count struct for profiling */
185 /* Local cycle count enum for profiling */
187 enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
190 /* Thread-local work struct, contains part of nbnxn_grid_t */
192 gmx_cache_protect_t cp0;
198 int sort_work_nalloc;
200 nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
202 int ndistc; /* Number of distance checks for flop counting */
204 nbnxn_cycle_t cc[enbsCCnr];
206 gmx_cache_protect_t cp1;
207 } nbnxn_search_work_t;
209 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
210 typedef struct nbnxn_search {
211 int ePBC; /* PBC type enum */
212 matrix box; /* The periodic unit-cell */
214 gmx_bool DomDec; /* Are we doing domain decomposition? */
215 ivec dd_dim; /* Are we doing DD in x,y,z? */
216 gmx_domdec_zones_t *zones; /* The domain decomposition zones */
218 int ngrid; /* The number of grids, equal to #DD-zones */
219 nbnxn_grid_t *grid; /* Array of grids, size ngrid */
220 int *cell; /* Actual allocated cell array for all grids */
221 int cell_nalloc; /* Allocation size of cell */
222 int *a; /* Atom index for grid, the inverse of cell */
223 int a_nalloc; /* Allocation size of a */
225 int natoms_local; /* The local atoms run from 0 to natoms_local */
226 int natoms_nonlocal; /* The non-local atoms run from natoms_local
227 * to natoms_nonlocal */
229 gmx_bool print_cycles;
231 nbnxn_cycle_t cc[enbsCCnr];
233 gmx_icell_set_x_t *icell_set_x; /* Function for setting i-coords */
235 int nthread_max; /* Maximum number of threads for pair-search */
236 nbnxn_search_work_t *work; /* Work array, size nthread_max */
240 static void nbs_cycle_start(nbnxn_cycle_t *cc)
242 cc->start = gmx_cycles_read();
245 static void nbs_cycle_stop(nbnxn_cycle_t *cc)
247 cc->c += gmx_cycles_read() - cc->start;