2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, 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.
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.
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.
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.
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.
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.
36 #ifndef _nbnxn_internal_h
37 #define _nbnxn_internal_h
42 #include "gromacs/domdec/domdec.h"
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/mdlib/nbnxn_pairlist.h"
45 #include "gromacs/simd/simd.h"
46 #include "gromacs/timing/cyclecounter.h"
47 #include "gromacs/utility/alignedallocator.h"
48 #include "gromacs/utility/arrayref.h"
49 #include "gromacs/utility/real.h"
51 struct gmx_domdec_zones_t;
54 /* The number of clusters in a pair-search cell, used for GPU */
55 static const int c_gpuNumClusterPerCellZ = 2;
56 static const int c_gpuNumClusterPerCellY = 2;
57 static const int c_gpuNumClusterPerCellX = 2;
58 static const int c_gpuNumClusterPerCell = c_gpuNumClusterPerCellZ*c_gpuNumClusterPerCellY*c_gpuNumClusterPerCellX;
61 /* Strides for x/f with xyz and xyzq coordinate (and charge) storage */
64 /* Size of packs of x, y or z with SIMD packed coords/forces */
65 static const int c_packX4 = 4;
66 static const int c_packX8 = 8;
67 /* Strides for a pack of 4 and 8 coordinates/forces */
68 #define STRIDE_P4 (DIM*c_packX4)
69 #define STRIDE_P8 (DIM*c_packX8)
71 /* Returns the index in a coordinate array corresponding to atom a */
72 template<int packSize> static inline int atom_to_x_index(int a)
74 return DIM*(a & ~(packSize - 1)) + (a & (packSize - 1));
79 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
80 #define NBNXN_MEM_ALIGN (GMX_SIMD_REAL_WIDTH*sizeof(real))
82 /* No alignment required, but set it so we can call the same routines */
83 #define NBNXN_MEM_ALIGN 32
87 /* Bounding box calculations are (currently) always in single precision, so
88 * we only need to check for single precision support here.
89 * This uses less (cache-)memory and SIMD is faster, at least on x86.
91 #if GMX_SIMD4_HAVE_FLOAT
92 # define NBNXN_SEARCH_BB_SIMD4 1
93 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
94 # define NBNXN_SEARCH_BB_MEM_ALIGN (GMX_SIMD4_WIDTH*sizeof(float))
96 # define NBNXN_SEARCH_BB_SIMD4 0
97 /* No alignment required, but set it so we can call the same routines */
98 # define NBNXN_SEARCH_BB_MEM_ALIGN 32
102 /* Pair search box lower and upper corner in x,y,z.
103 * Store this in 4 iso 3 reals, which is useful with 4-wide SIMD.
104 * To avoid complicating the code we also use 4 without 4-wide SIMD.
107 /* Pair search box lower and upper bound in z only. */
109 /* Pair search box lower and upper corner x,y,z indices, entry 3 is unused */
115 #if NBNXN_SEARCH_BB_SIMD4
116 /* Always use 4-wide SIMD for bounding box calculations */
119 /* Single precision BBs + coordinates, we can also load coordinates with SIMD */
120 # define NBNXN_SEARCH_SIMD4_FLOAT_X_BB 1
122 # define NBNXN_SEARCH_SIMD4_FLOAT_X_BB 0
125 /* The packed bounding box coordinate stride is always set to 4.
126 * With AVX we could use 8, but that turns out not to be faster.
128 # define STRIDE_PBB GMX_SIMD4_WIDTH
129 # define STRIDE_PBB_2LOG 2
131 /* Store bounding boxes corners as quadruplets: xxxxyyyyzzzz */
132 # define NBNXN_BBXXXX 1
133 /* Size of bounding box corners quadruplet */
134 # define NNBSBB_XXXX (NNBSBB_D*DIM*STRIDE_PBB)
136 #else /* NBNXN_SEARCH_BB_SIMD4 */
138 # define NBNXN_SEARCH_SIMD4_FLOAT_X_BB 0
139 # define NBNXN_BBXXXX 0
141 #endif /* NBNXN_SEARCH_BB_SIMD4 */
144 /* Bounding box for a nbnxn atom cluster */
146 float lower[NNBSBB_C];
147 float upper[NNBSBB_C];
151 /* A pair-search grid struct for one domain decomposition zone */
154 rvec c0; /* The lower corner of the (local) grid */
155 rvec c1; /* The upper corner of the (local) grid */
156 rvec size; /* c1 - c0 */
157 real atom_density; /* The atom number density for the local grid */
159 gmx_bool bSimple; /* Is this grid simple or super/sub */
160 int na_c; /* Number of atoms per cluster */
161 int na_cj; /* Number of atoms for list j-clusters */
162 int na_sc; /* Number of atoms per super-cluster */
163 int na_c_2log; /* 2log of na_c */
165 int numCells[DIM - 1]; /* Number of cells along x/y */
166 int nc; /* Total number of cells */
168 real cellSize[DIM - 1]; /* size of a cell */
169 real invCellSize[DIM - 1]; /* 1/cellSize */
171 int cell0; /* Index in nbs->cell corresponding to cell 0 */
174 std::vector<int> cxy_na; /* The number of atoms for each column in x,y */
175 std::vector<int> cxy_ind; /* Grid (super)cell index, offset from cell0 */
177 std::vector<int> nsubc; /* The number of sub cells for each super cell */
180 std::vector<float> bbcz; /* Bounding boxes in z for the cells */
181 std::vector < nbnxn_bb_t, gmx::AlignedAllocator < nbnxn_bb_t>> bb; /* 3D bounding boxes for the sub cells */
182 std::vector < nbnxn_bb_t, gmx::AlignedAllocator < nbnxn_bb_t>> bbjStorage; /* 3D j-bounding boxes for the case where
183 * the i- and j-cluster sizes are different */
184 gmx::ArrayRef<nbnxn_bb_t> bbj; /* 3D j-bounding boxes */
185 std::vector < float, gmx::AlignedAllocator < float>> pbb; /* 3D b. boxes in xxxx format per super cell */
187 /* Bit-flag information */
188 std::vector<int> flags; /* Flags for properties of clusters in each cell */
189 std::vector<unsigned int> fep; /* FEP signal bits for atoms in each cluster */
191 int nsubc_tot; /* Total number of subcell, used for printing */
194 /* Working data for the actual i-supercell during pair search */
195 typedef struct nbnxn_list_work {
196 gmx_cache_protect_t cp0; /* Protect cache between threads */
198 nbnxn_bb_t *bb_ci; /* The bounding boxes, pbc shifted, for each cluster */
199 float *pbb_ci; /* As bb_ci, but in xxxx packed format */
200 real *x_ci; /* The coordinates, pbc shifted, for each atom */
201 real *x_ci_simd; /* aligned pointer to 4*DIM*GMX_SIMD_REAL_WIDTH floats */
202 int cj_ind; /* The current cj_ind index for the current list */
203 int cj4_init; /* The first unitialized cj4 block */
205 float *d2; /* Bounding box distance work array */
207 nbnxn_cj_t *cj; /* The j-cell list */
208 int cj_nalloc; /* Allocation size of cj */
210 int ncj_noq; /* Nr. of cluster pairs without Coul for flop count */
211 int ncj_hlj; /* Nr. of cluster pairs with 1/2 LJ for flop count */
213 int *sort; /* Sort index */
214 int sort_nalloc; /* Allocation size of sort */
216 nbnxn_sci_t *sci_sort; /* Second sci array, for sorting */
217 int sci_sort_nalloc; /* Allocation size of sci_sort */
219 gmx_cache_protect_t cp1; /* Protect cache between threads */
222 /* Function type for setting the i-atom coordinate working data */
224 gmx_icell_set_x_t (int ci,
225 real shx, real shy, real shz,
226 int stride, const real *x,
227 nbnxn_list_work_t *work);
229 /* Local cycle count struct for profiling */
236 /* Local cycle count enum for profiling */
238 enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
241 /* Thread-local work struct, contains part of nbnxn_grid_t */
242 struct nbnxn_search_work_t
245 nbnxn_search_work_t();
248 ~nbnxn_search_work_t();
250 gmx_cache_protect_t cp0; /* Buffer to avoid cache polution */
252 std::vector<int> cxy_na; /* Grid column atom counts temporary buffer */
254 std::vector<int> sortBuffer; /* Temporary buffer for sorting atoms within a grid column */
256 nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
258 int ndistc; /* Number of distance checks for flop counting */
261 std::unique_ptr<t_nblist> nbl_fep; /* Temporary FEP list for load balancing */
263 nbnxn_cycle_t cc[enbsCCnr]; /* Counters for thread-local cycles */
265 gmx_cache_protect_t cp1; /* Buffer to avoid cache polution */
268 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
271 /* \brief Constructor
273 * \param[in] n_dd_cells The number of domain decomposition cells per dimension, without DD nullptr should be passed
274 * \param[in] zones The domain decomposition zone setup, without DD nullptr should be passed
275 * \param[in] bFEP Tells whether non-bonded interactions are perturbed
276 * \param[in] nthread_max The maximum number of threads used in the search
279 nbnxn_search(const ivec *n_dd_cells,
280 const gmx_domdec_zones_t *zones,
284 gmx_bool bFEP; /* Do we have perturbed atoms? */
285 int ePBC; /* PBC type enum */
286 matrix box; /* The periodic unit-cell */
288 gmx_bool DomDec; /* Are we doing domain decomposition? */
289 ivec dd_dim; /* Are we doing DD in x,y,z? */
290 const gmx_domdec_zones_t *zones; /* The domain decomposition zones */
292 std::vector<nbnxn_grid_t> grid; /* Array of grids, size ngrid */
293 std::vector<int> cell; /* Actual allocated cell array for all grids */
294 std::vector<int> a; /* Atom index for grid, the inverse of cell */
296 int natoms_local; /* The local atoms run from 0 to natoms_local */
297 int natoms_nonlocal; /* The non-local atoms run from natoms_local
298 * to natoms_nonlocal */
300 gmx_bool print_cycles;
302 nbnxn_cycle_t cc[enbsCCnr];
304 gmx_icell_set_x_t *icell_set_x; /* Function for setting i-coords */
306 /* Thread-local work data */
307 mutable std::vector<nbnxn_search_work_t> work; /* Work array, one entry for each thread */
311 static inline void nbs_cycle_start(nbnxn_cycle_t *cc)
313 cc->start = gmx_cycles_read();
316 static inline void nbs_cycle_stop(nbnxn_cycle_t *cc)
318 cc->c += gmx_cycles_read() - cc->start;