2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, 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 /* Strides for x/f with xyz and xyzq coordinate (and charge) storage */
57 /* Size of packs of x, y or z with SIMD packed coords/forces */
58 static const int c_packX4 = 4;
59 static const int c_packX8 = 8;
60 /* Strides for a pack of 4 and 8 coordinates/forces */
61 #define STRIDE_P4 (DIM*c_packX4)
62 #define STRIDE_P8 (DIM*c_packX8)
64 /* Returns the index in a coordinate array corresponding to atom a */
65 template<int packSize> static inline int atom_to_x_index(int a)
67 return DIM*(a & ~(packSize - 1)) + (a & (packSize - 1));
72 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
73 #define NBNXN_MEM_ALIGN (GMX_SIMD_REAL_WIDTH*sizeof(real))
75 /* No alignment required, but set it so we can call the same routines */
76 #define NBNXN_MEM_ALIGN 32
80 /* Bounding box calculations are (currently) always in single precision, so
81 * we only need to check for single precision support here.
82 * This uses less (cache-)memory and SIMD is faster, at least on x86.
84 #if GMX_SIMD4_HAVE_FLOAT
85 # define NBNXN_SEARCH_BB_SIMD4 1
86 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
87 # define NBNXN_SEARCH_BB_MEM_ALIGN (GMX_SIMD4_WIDTH*sizeof(float))
89 # define NBNXN_SEARCH_BB_SIMD4 0
90 /* No alignment required, but set it so we can call the same routines */
91 # define NBNXN_SEARCH_BB_MEM_ALIGN 32
95 /* Pair search box lower and upper corner in x,y,z.
96 * Store this in 4 iso 3 reals, which is useful with 4-wide SIMD.
97 * To avoid complicating the code we also use 4 without 4-wide SIMD.
100 /* Pair search box lower and upper bound in z only. */
102 /* Pair search box lower and upper corner x,y,z indices, entry 3 is unused */
108 #if NBNXN_SEARCH_BB_SIMD4
109 /* Always use 4-wide SIMD for bounding box calculations */
112 /* Single precision BBs + coordinates, we can also load coordinates with SIMD */
113 # define NBNXN_SEARCH_SIMD4_FLOAT_X_BB 1
115 # define NBNXN_SEARCH_SIMD4_FLOAT_X_BB 0
118 /* The packed bounding box coordinate stride is always set to 4.
119 * With AVX we could use 8, but that turns out not to be faster.
121 # define STRIDE_PBB GMX_SIMD4_WIDTH
122 # define STRIDE_PBB_2LOG 2
124 /* Store bounding boxes corners as quadruplets: xxxxyyyyzzzz */
125 # define NBNXN_BBXXXX 1
126 /* Size of bounding box corners quadruplet */
127 # define NNBSBB_XXXX (NNBSBB_D*DIM*STRIDE_PBB)
129 #else /* NBNXN_SEARCH_BB_SIMD4 */
131 # define NBNXN_SEARCH_SIMD4_FLOAT_X_BB 0
132 # define NBNXN_BBXXXX 0
134 #endif /* NBNXN_SEARCH_BB_SIMD4 */
138 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
141 /* Bounding box for a nbnxn atom cluster */
143 float lower[NNBSBB_C];
144 float upper[NNBSBB_C];
148 /* A pair-search grid struct for one domain decomposition zone
150 * Note that when atom groups, instead of individual atoms, are assigned
151 * to grid cells, individual atoms can be geometrically outside the cell
152 * and grid that they have been assigned to (as determined by the center
153 * or geometry of the atom group they belong to).
157 rvec c0; /* The lower corner of the (local) grid */
158 rvec c1; /* The upper corner of the (local) grid */
159 rvec size; /* c1 - c0 */
160 real atom_density; /* The atom number density for the local grid */
161 real maxAtomGroupRadius; /* The maximum distance an atom can be outside
162 * of a cell and outside of the grid
165 gmx_bool bSimple; /* Is this grid simple or super/sub */
166 int na_c; /* Number of atoms per cluster */
167 int na_cj; /* Number of atoms for list j-clusters */
168 int na_sc; /* Number of atoms per super-cluster */
169 int na_c_2log; /* 2log of na_c */
171 int numCells[DIM - 1]; /* Number of cells along x/y */
172 int nc; /* Total number of cells */
174 real cellSize[DIM - 1]; /* size of a cell */
175 real invCellSize[DIM - 1]; /* 1/cellSize */
177 int cell0; /* Index in nbs->cell corresponding to cell 0 */
180 std::vector<int> cxy_na; /* The number of atoms for each column in x,y */
181 std::vector<int> cxy_ind; /* Grid (super)cell index, offset from cell0 */
183 std::vector<int> nsubc; /* The number of sub cells for each super cell */
186 std::vector<float> bbcz; /* Bounding boxes in z for the cells */
187 std::vector < nbnxn_bb_t, gmx::AlignedAllocator < nbnxn_bb_t>> bb; /* 3D bounding boxes for the sub cells */
188 std::vector < nbnxn_bb_t, gmx::AlignedAllocator < nbnxn_bb_t>> bbjStorage; /* 3D j-bounding boxes for the case where
189 * the i- and j-cluster sizes are different */
190 gmx::ArrayRef<nbnxn_bb_t> bbj; /* 3D j-bounding boxes */
191 std::vector < float, gmx::AlignedAllocator < float>> pbb; /* 3D b. boxes in xxxx format per super cell */
193 /* Bit-flag information */
194 std::vector<int> flags; /* Flags for properties of clusters in each cell */
195 std::vector<unsigned int> fep; /* FEP signal bits for atoms in each cluster */
198 int nsubc_tot; /* Total number of subcell, used for printing */
201 /* Working data for the actual i-supercell during pair search */
202 struct NbnxnPairlistCpuWork
204 // Struct for storing coordinats and bounding box for an i-entry during search
209 x(c_nbnxnCpuIClusterSize*DIM),
210 xSimd(c_nbnxnCpuIClusterSize*DIM*GMX_REAL_MAX_SIMD_WIDTH)
214 // The bounding boxes, pbc shifted, for each cluster
215 AlignedVector<nbnxn_bb_t> bb;
216 // The coordinates, pbc shifted, for each atom
218 // Aligned list for storing 4*DIM*GMX_SIMD_REAL_WIDTH reals
219 AlignedVector<real> xSimd;
222 // Protect data from cache pollution between threads
223 gmx_cache_protect_t cp0;
225 // Work data for generating an IEntry in the pairlist
226 IClusterData iClusterData;
227 // The current cj_ind index for the current list
229 // Temporary j-cluster list, used for sorting on exclusions
230 std::vector<nbnxn_cj_t> cj;
232 // Nr. of cluster pairs without Coulomb for flop counting
234 // Nr. of cluster pairs with 1/2 LJ for flop count
237 // Protect data from cache pollution between threads
238 gmx_cache_protect_t cp1;
241 /* Working data for the actual i-supercell during pair search */
242 struct NbnxnPairlistGpuWork
244 struct ISuperClusterData
246 ISuperClusterData() :
247 bb(c_gpuNumClusterPerCell),
248 #if NBNXN_SEARCH_BB_SIMD4
249 bbPacked(c_gpuNumClusterPerCell/STRIDE_PBB*NNBSBB_XXXX),
251 x(c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize*DIM),
252 xSimd(c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize*DIM)
256 // The bounding boxes, pbc shifted, for each cluster
257 AlignedVector<nbnxn_bb_t> bb;
258 // As bb, but in packed xxxx format
259 AlignedVector<float> bbPacked;
260 // The coordinates, pbc shifted, for each atom
261 AlignedVector<real> x;
262 // Aligned coordinate list used for 4*DIM*GMX_SIMD_REAL_WIDTH floats
263 AlignedVector<real> xSimd;
266 NbnxnPairlistGpuWork() :
267 distanceBuffer(c_gpuNumClusterPerCell),
268 sci_sort({}, {gmx::PinningPolicy::PinnedIfSupported})
272 // Protect data from cache pollution between threads
273 gmx_cache_protect_t cp0;
275 // Work data for generating an i-entry in the pairlist
276 ISuperClusterData iSuperClusterData;
277 // The current j-cluster index for the current list
279 // Bounding box distance work array
280 AlignedVector<float> distanceBuffer;
282 // Buffer for sorting list entries
283 std::vector<int> sortBuffer;
285 // Second sci array, for sorting
286 gmx::HostVector<nbnxn_sci_t> sci_sort;
288 // Protect data from cache pollution between threads
289 gmx_cache_protect_t cp1;
292 /* Local cycle count struct for profiling */
299 /* Local cycle count enum for profiling */
301 enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
304 /* Thread-local work struct, contains part of nbnxn_grid_t */
305 struct nbnxn_search_work_t
307 nbnxn_search_work_t();
309 ~nbnxn_search_work_t();
311 gmx_cache_protect_t cp0; /* Buffer to avoid cache polution */
313 std::vector<int> cxy_na; /* Grid column atom counts temporary buffer */
315 std::vector<int> sortBuffer; /* Temporary buffer for sorting atoms within a grid column */
317 nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
319 int ndistc; /* Number of distance checks for flop counting */
322 std::unique_ptr<t_nblist> nbl_fep; /* Temporary FEP list for load balancing */
324 nbnxn_cycle_t cc[enbsCCnr]; /* Counters for thread-local cycles */
326 gmx_cache_protect_t cp1; /* Buffer to avoid cache polution */
329 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
332 /* \brief Constructor
334 * \param[in] n_dd_cells The number of domain decomposition cells per dimension, without DD nullptr should be passed
335 * \param[in] zones The domain decomposition zone setup, without DD nullptr should be passed
336 * \param[in] bFEP Tells whether non-bonded interactions are perturbed
337 * \param[in] nthread_max The maximum number of threads used in the search
340 nbnxn_search(const ivec *n_dd_cells,
341 const gmx_domdec_zones_t *zones,
345 gmx_bool bFEP; /* Do we have perturbed atoms? */
346 int ePBC; /* PBC type enum */
347 matrix box; /* The periodic unit-cell */
349 gmx_bool DomDec; /* Are we doing domain decomposition? */
350 ivec dd_dim; /* Are we doing DD in x,y,z? */
351 const gmx_domdec_zones_t *zones; /* The domain decomposition zones */
353 std::vector<nbnxn_grid_t> grid; /* Array of grids, size ngrid */
354 std::vector<int> cell; /* Actual allocated cell array for all grids */
355 std::vector<int> a; /* Atom index for grid, the inverse of cell */
357 int natoms_local; /* The local atoms run from 0 to natoms_local */
358 int natoms_nonlocal; /* The non-local atoms run from natoms_local
359 * to natoms_nonlocal */
361 gmx_bool print_cycles;
363 nbnxn_cycle_t cc[enbsCCnr];
365 /* Thread-local work data */
366 mutable std::vector<nbnxn_search_work_t> work; /* Work array, one entry for each thread */
370 static inline void nbs_cycle_start(nbnxn_cycle_t *cc)
372 cc->start = gmx_cycles_read();
375 static inline void nbs_cycle_stop(nbnxn_cycle_t *cc)
377 cc->c += gmx_cycles_read() - cc->start;