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
41 #include "gromacs/domdec/domdec.h"
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/mdlib/nbnxn_pairlist.h"
44 #include "gromacs/simd/simd.h"
45 #include "gromacs/timing/cyclecounter.h"
46 #include "gromacs/utility/alignedallocator.h"
47 #include "gromacs/utility/arrayref.h"
48 #include "gromacs/utility/real.h"
50 using namespace gmx
; // TODO: Remove when this file is moved into gmx namespace
52 struct gmx_domdec_zones_t
;
55 /* The number of clusters in a pair-search cell, used for GPU */
56 static const int c_gpuNumClusterPerCellZ
= 2;
57 static const int c_gpuNumClusterPerCellY
= 2;
58 static const int c_gpuNumClusterPerCellX
= 2;
59 static const int c_gpuNumClusterPerCell
= c_gpuNumClusterPerCellZ
*c_gpuNumClusterPerCellY
*c_gpuNumClusterPerCellX
;
62 /* Strides for x/f with xyz and xyzq coordinate (and charge) storage */
65 /* Size of packs of x, y or z with SIMD packed coords/forces */
66 static const int c_packX4
= 4;
67 static const int c_packX8
= 8;
68 /* Strides for a pack of 4 and 8 coordinates/forces */
69 #define STRIDE_P4 (DIM*c_packX4)
70 #define STRIDE_P8 (DIM*c_packX8)
72 /* Returns the index in a coordinate array corresponding to atom a */
73 template<int packSize
> static inline int atom_to_x_index(int a
)
75 return DIM
*(a
& ~(packSize
- 1)) + (a
& (packSize
- 1));
80 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
81 #define NBNXN_MEM_ALIGN (GMX_SIMD_REAL_WIDTH*sizeof(real))
83 /* No alignment required, but set it so we can call the same routines */
84 #define NBNXN_MEM_ALIGN 32
88 /* Bounding box calculations are (currently) always in single precision, so
89 * we only need to check for single precision support here.
90 * This uses less (cache-)memory and SIMD is faster, at least on x86.
92 #if GMX_SIMD4_HAVE_FLOAT
93 # define NBNXN_SEARCH_BB_SIMD4 1
94 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
95 # define NBNXN_SEARCH_BB_MEM_ALIGN (GMX_SIMD4_WIDTH*sizeof(float))
97 # define NBNXN_SEARCH_BB_SIMD4 0
98 /* No alignment required, but set it so we can call the same routines */
99 # define NBNXN_SEARCH_BB_MEM_ALIGN 32
103 /* Pair search box lower and upper corner in x,y,z.
104 * Store this in 4 iso 3 reals, which is useful with 4-wide SIMD.
105 * To avoid complicating the code we also use 4 without 4-wide SIMD.
108 /* Pair search box lower and upper bound in z only. */
110 /* Pair search box lower and upper corner x,y,z indices, entry 3 is unused */
116 #if NBNXN_SEARCH_BB_SIMD4
117 /* Always use 4-wide SIMD for bounding box calculations */
120 /* Single precision BBs + coordinates, we can also load coordinates with SIMD */
121 # define NBNXN_SEARCH_SIMD4_FLOAT_X_BB 1
123 # define NBNXN_SEARCH_SIMD4_FLOAT_X_BB 0
126 /* The packed bounding box coordinate stride is always set to 4.
127 * With AVX we could use 8, but that turns out not to be faster.
129 # define STRIDE_PBB GMX_SIMD4_WIDTH
130 # define STRIDE_PBB_2LOG 2
132 /* Store bounding boxes corners as quadruplets: xxxxyyyyzzzz */
133 # define NBNXN_BBXXXX 1
134 /* Size of bounding box corners quadruplet */
135 # define NNBSBB_XXXX (NNBSBB_D*DIM*STRIDE_PBB)
137 #else /* NBNXN_SEARCH_BB_SIMD4 */
139 # define NBNXN_SEARCH_SIMD4_FLOAT_X_BB 0
140 # define NBNXN_BBXXXX 0
142 #endif /* NBNXN_SEARCH_BB_SIMD4 */
145 /* Bounding box for a nbnxn atom cluster */
147 float lower
[NNBSBB_C
];
148 float upper
[NNBSBB_C
];
152 /* A pair-search grid struct for one domain decomposition zone */
155 rvec c0
; /* The lower corner of the (local) grid */
156 rvec c1
; /* The upper corner of the (local) grid */
157 rvec size
; /* c1 - c0 */
158 real atom_density
; /* The atom number density for the local grid */
160 gmx_bool bSimple
; /* Is this grid simple or super/sub */
161 int na_c
; /* Number of atoms per cluster */
162 int na_cj
; /* Number of atoms for list j-clusters */
163 int na_sc
; /* Number of atoms per super-cluster */
164 int na_c_2log
; /* 2log of na_c */
166 int ncx
; /* Number of (super-)cells along x */
167 int ncy
; /* Number of (super-)cells along y */
168 int nc
; /* Total number of (super-)cells */
170 real sx
; /* x-size of a (super-)cell */
171 real sy
; /* y-size of a (super-)cell */
172 real inv_sx
; /* 1/sx */
173 real inv_sy
; /* 1/sy */
175 int cell0
; /* Index in nbs->cell corresponding to cell 0 */
178 std::vector
<int> cxy_na
; /* The number of atoms for each column in x,y */
179 std::vector
<int> cxy_ind
; /* Grid (super)cell index, offset from cell0 */
181 std::vector
<int> nsubc
; /* The number of sub cells for each super cell */
184 std::vector
<float> bbcz
; /* Bounding boxes in z for the cells */
185 std::vector
< nbnxn_bb_t
, AlignedAllocator
< nbnxn_bb_t
>> bb
; /* 3D bounding boxes for the sub cells */
186 std::vector
< nbnxn_bb_t
, AlignedAllocator
< nbnxn_bb_t
>> bbjStorage
; /* 3D j-bounding boxes for the case where
187 * the i- and j-cluster sizes are different */
188 gmx::ArrayRef
<nbnxn_bb_t
> bbj
; /* 3D j-bounding boxes */
189 std::vector
< float, AlignedAllocator
< float>> pbb
; /* 3D b. boxes in xxxx format per super cell */
191 /* Bit-flag information */
192 std::vector
<int> flags
; /* Flags for properties of clusters in each cell */
193 std::vector
<unsigned int> fep
; /* FEP signal bits for atoms in each cluster */
195 int nsubc_tot
; /* Total number of subcell, used for printing */
198 /* Working data for the actual i-supercell during pair search */
199 typedef struct nbnxn_list_work
{
200 gmx_cache_protect_t cp0
; /* Protect cache between threads */
202 nbnxn_bb_t
*bb_ci
; /* The bounding boxes, pbc shifted, for each cluster */
203 float *pbb_ci
; /* As bb_ci, but in xxxx packed format */
204 real
*x_ci
; /* The coordinates, pbc shifted, for each atom */
205 real
*x_ci_simd
; /* aligned pointer to 4*DIM*GMX_SIMD_REAL_WIDTH floats */
206 int cj_ind
; /* The current cj_ind index for the current list */
207 int cj4_init
; /* The first unitialized cj4 block */
209 float *d2
; /* Bounding box distance work array */
211 nbnxn_cj_t
*cj
; /* The j-cell list */
212 int cj_nalloc
; /* Allocation size of cj */
214 int ncj_noq
; /* Nr. of cluster pairs without Coul for flop count */
215 int ncj_hlj
; /* Nr. of cluster pairs with 1/2 LJ for flop count */
217 int *sort
; /* Sort index */
218 int sort_nalloc
; /* Allocation size of sort */
220 nbnxn_sci_t
*sci_sort
; /* Second sci array, for sorting */
221 int sci_sort_nalloc
; /* Allocation size of sci_sort */
223 gmx_cache_protect_t cp1
; /* Protect cache between threads */
226 /* Function type for setting the i-atom coordinate working data */
228 gmx_icell_set_x_t (int ci
,
229 real shx
, real shy
, real shz
,
230 int stride
, const real
*x
,
231 nbnxn_list_work_t
*work
);
233 /* Local cycle count struct for profiling */
240 /* Local cycle count enum for profiling */
242 enbsCCgrid
, enbsCCsearch
, enbsCCcombine
, enbsCCreducef
, enbsCCnr
245 /* Thread-local work struct, contains part of nbnxn_grid_t */
247 gmx_cache_protect_t cp0
;
253 int sort_work_nalloc
;
255 nbnxn_buffer_flags_t buffer_flags
; /* Flags for force buffer access */
257 int ndistc
; /* Number of distance checks for flop counting */
259 t_nblist
*nbl_fep
; /* Temporary FEP list for load balancing */
261 nbnxn_cycle_t cc
[enbsCCnr
];
263 gmx_cache_protect_t cp1
;
264 } nbnxn_search_work_t
;
266 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
267 typedef struct nbnxn_search
{
268 gmx_bool bFEP
; /* Do we have perturbed atoms? */
269 int ePBC
; /* PBC type enum */
270 matrix box
; /* The periodic unit-cell */
272 gmx_bool DomDec
; /* Are we doing domain decomposition? */
273 ivec dd_dim
; /* Are we doing DD in x,y,z? */
274 struct gmx_domdec_zones_t
*zones
; /* The domain decomposition zones */
276 std::vector
<nbnxn_grid_t
> grid
; /* Array of grids, size ngrid */
277 int *cell
; /* Actual allocated cell array for all grids */
278 int cell_nalloc
; /* Allocation size of cell */
279 int *a
; /* Atom index for grid, the inverse of cell */
280 int a_nalloc
; /* Allocation size of a */
282 int natoms_local
; /* The local atoms run from 0 to natoms_local */
283 int natoms_nonlocal
; /* The non-local atoms run from natoms_local
284 * to natoms_nonlocal */
286 gmx_bool print_cycles
;
288 nbnxn_cycle_t cc
[enbsCCnr
];
290 gmx_icell_set_x_t
*icell_set_x
; /* Function for setting i-coords */
292 int nthread_max
; /* Maximum number of threads for pair-search */
293 nbnxn_search_work_t
*work
; /* Work array, size nthread_max */
297 static inline void nbs_cycle_start(nbnxn_cycle_t
*cc
)
299 cc
->start
= gmx_cycles_read();
302 static inline void nbs_cycle_stop(nbnxn_cycle_t
*cc
)
304 cc
->c
+= gmx_cycles_read() - cc
->start
;