2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017, 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_pairlist_h
37 #define _nbnxn_pairlist_h
41 #include "gromacs/math/vectypes.h"
42 #include "gromacs/mdtypes/nblist.h"
43 #include "gromacs/utility/basedefinitions.h"
44 #include "gromacs/utility/bitmask.h"
45 #include "gromacs/utility/real.h"
49 /* With GPU kernels the i and j cluster size is 8 atoms */
50 static const int c_nbnxnGpuClusterSize = 8;
52 /* The number of clusters in a super-cluster, used for GPU */
53 static const int c_nbnxnGpuNumClusterPerSupercluster = 8;
55 /* With GPU kernels we group cluster pairs in 4 to optimize memory usage
56 * of integers containing 32 bits.
58 static const int c_nbnxnGpuJgroupSize = 32/c_nbnxnGpuNumClusterPerSupercluster;
60 /* In CUDA the number of threads in a warp is 32 and we have cluster pairs
61 * of 8*8=64 atoms, so it's convenient to store data for cluster pair halves.
63 static const int c_nbnxnGpuClusterpairSplit = 2;
65 /* The fixed size of the exclusion mask array for a half cluster pair */
66 static const int c_nbnxnGpuExclSize = c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
68 /* A buffer data structure of 64 bytes
69 * to be placed at the beginning and end of structs
70 * to avoid cache invalidation of the real contents
71 * of the struct by writes to neighboring memory.
75 } gmx_cache_protect_t;
77 /* Abstract type for pair searching data */
78 typedef struct nbnxn_search * nbnxn_search_t;
80 /* Function that should return a pointer *ptr to memory
82 * Error handling should be done within this function.
84 typedef void nbnxn_alloc_t (void **ptr, size_t nbytes);
86 /* Function that should free the memory pointed to by *ptr.
87 * NULL should not be passed to this function.
89 typedef void nbnxn_free_t (void *ptr);
91 /* This is the actual cluster-pair list j-entry.
92 * cj is the j-cluster.
93 * The interaction bits in excl are indexed i-major, j-minor.
94 * The cj entries are sorted such that ones with exclusions come first.
95 * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
96 * is found, all subsequent j-entries in the i-entry also have full masks.
99 int cj; /* The j-cluster */
100 unsigned int excl; /* The exclusion (interaction) bits */
103 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
104 * The upper bits contain information for non-bonded kernel optimization.
105 * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
106 * But three flags can be used to skip interactions, currently only for subc=0
107 * !(shift & NBNXN_CI_DO_LJ(subc)) => we can skip LJ for all pairs
108 * shift & NBNXN_CI_HALF_LJ(subc) => we can skip LJ for the second half of i
109 * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
111 #define NBNXN_CI_SHIFT 127
112 #define NBNXN_CI_DO_LJ(subc) (1<<(7+3*(subc)))
113 #define NBNXN_CI_HALF_LJ(subc) (1<<(8+3*(subc)))
114 #define NBNXN_CI_DO_COUL(subc) (1<<(9+3*(subc)))
116 /* Simple pair-list i-unit */
118 int ci; /* i-cluster */
119 int shift; /* Shift vector index plus possible flags, see above */
120 int cj_ind_start; /* Start index into cj */
121 int cj_ind_end; /* End index into cj */
124 /* Grouped pair-list i-unit */
126 int sci; /* i-super-cluster */
127 int shift; /* Shift vector index plus possible flags */
128 int cj4_ind_start; /* Start index into cj4 */
129 int cj4_ind_end; /* End index into cj4 */
133 unsigned int imask; /* The i-cluster interactions mask for 1 warp */
134 int excl_ind; /* Index into the exclusion array for 1 warp */
138 int cj[c_nbnxnGpuJgroupSize]; /* The 4 j-clusters */
139 nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data for 2 warps */
143 unsigned int pair[c_nbnxnGpuExclSize]; /* Topology exclusion interaction bits for one warp,
144 * each unsigned has bitS for 4*8 i clusters
148 typedef struct nbnxn_pairlist_t {
149 gmx_cache_protect_t cp0;
151 nbnxn_alloc_t *alloc;
154 gmx_bool bSimple; /* Simple list has na_sc=na_s and uses cj *
155 * Complex list uses cj4 */
157 int na_ci; /* The number of atoms per i-cluster */
158 int na_cj; /* The number of atoms per j-cluster */
159 int na_sc; /* The number of atoms per super cluster */
160 real rlist; /* The radius for constructing the list */
161 int nci; /* The number of i-clusters in the list */
162 nbnxn_ci_t *ci; /* The i-cluster list, size nci */
163 int ci_nalloc; /* The allocation size of ci */
164 int nsci; /* The number of i-super-clusters in the list */
165 nbnxn_sci_t *sci; /* The i-super-cluster list */
166 int sci_nalloc; /* The allocation size of sci */
168 int ncj; /* The number of j-clusters in the list */
169 nbnxn_cj_t *cj; /* The j-cluster list, size ncj */
170 int cj_nalloc; /* The allocation size of cj */
171 int ncjInUse; /* The number of j-clusters that are used by ci entries in this list, will be <= ncj */
173 int ncj4; /* The total number of 4*j clusters */
174 nbnxn_cj4_t *cj4; /* The 4*j cluster list, size ncj4 */
175 int cj4_nalloc; /* The allocation size of cj4 */
176 int nexcl; /* The count for excl */
177 nbnxn_excl_t *excl; /* Atom interaction bits (non-exclusions) */
178 int excl_nalloc; /* The allocation size for excl */
179 int nci_tot; /* The total number of i clusters */
181 struct nbnxn_list_work *work;
183 gmx_cache_protect_t cp1;
187 int nnbl; /* number of lists */
188 nbnxn_pairlist_t **nbl; /* lists */
189 nbnxn_pairlist_t **nbl_work; /* work space for rebalancing lists */
190 gmx_bool bCombined; /* TRUE if lists get combined into one (the 1st) */
191 gmx_bool bSimple; /* TRUE if the list of of type "simple"
192 (na_sc=na_s, no super-clusters used) */
193 int natpair_ljq; /* Total number of atom pairs for LJ+Q kernel */
194 int natpair_lj; /* Total number of atom pairs for LJ kernel */
195 int natpair_q; /* Total number of atom pairs for Q kernel */
197 } nbnxn_pairlist_set_t;
200 nbatXYZ, nbatXYZQ, nbatX4, nbatX8
204 real *f; /* f, size natoms*fstride */
205 real *fshift; /* Shift force array, size SHIFTS*DIM */
206 int nV; /* The size of *Vvdw and *Vc */
207 real *Vvdw; /* Temporary Van der Waals group energy storage */
208 real *Vc; /* Temporary Coulomb group energy storage */
209 int nVS; /* The size of *VSvdw and *VSc */
210 real *VSvdw; /* Temporary SIMD Van der Waals group energy storage */
211 real *VSc; /* Temporary SIMD Coulomb group energy storage */
212 } nbnxn_atomdata_output_t;
214 /* Block size in atoms for the non-bonded thread force-buffer reduction,
215 * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
216 * Should be small to reduce the reduction and zeroing cost,
217 * but too small will result in overhead.
218 * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
221 #define NBNXN_BUFFERFLAG_SIZE 8
223 #define NBNXN_BUFFERFLAG_SIZE 16
226 /* We store the reduction flags as gmx_bitmask_t.
227 * This limits the number of flags to BITMASK_SIZE.
229 #define NBNXN_BUFFERFLAG_MAX_THREADS (BITMASK_SIZE)
231 /* Flags for telling if threads write to force output buffers */
233 int nflag; /* The number of flag blocks */
234 gmx_bitmask_t *flag; /* Bit i is set when thread i writes to a cell-block */
235 int flag_nalloc; /* Allocation size of cxy_flag */
236 } nbnxn_buffer_flags_t;
238 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
240 ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
243 typedef struct nbnxn_atomdata_t {
244 nbnxn_alloc_t *alloc;
246 int ntype; /* The number of different atom types */
247 real *nbfp; /* Lennard-Jones 6*C6 and 12*C12 params, size ntype^2*2 */
248 int comb_rule; /* Combination rule, see enum above */
249 real *nbfp_comb; /* LJ parameter per atom type, size ntype*2 */
250 real *nbfp_aligned; /* As nbfp, but with an alignment (stride) suitable
251 * for the present SIMD architectures
253 int natoms; /* Number of atoms */
254 int natoms_local; /* Number of local atoms */
255 int *type; /* Atom types */
256 real *lj_comb; /* LJ parameters per atom for combining for pairs */
257 int XFormat; /* The format of x (and q), enum */
258 int FFormat; /* The format of f, enum */
259 real *q; /* Charges, can be NULL if incorporated in x */
260 int na_c; /* The number of atoms per cluster */
261 int nenergrp; /* The number of energy groups */
262 int neg_2log; /* Log2 of nenergrp */
263 int *energrp; /* The energy groups per cluster, can be NULL */
264 gmx_bool bDynamicBox; /* Do we need to update shift_vec every step? */
265 rvec *shift_vec; /* Shift vectors, copied from t_forcerec */
266 int xstride; /* stride for a coordinate in x (usually 3 or 4) */
267 int fstride; /* stride for a coordinate in f (usually 3 or 4) */
268 real *x; /* x and possibly q, size natoms*xstride */
270 /* j-atom minus i-atom index for generating self and Newton exclusions
271 * cluster-cluster pairs of the diagonal, for 4xn and 2xnn kernels.
273 real *simd_4xn_diagonal_j_minus_i;
274 real *simd_2xnn_diagonal_j_minus_i;
275 /* Filters for topology exclusion masks for the SIMD kernels. */
276 gmx_uint32_t *simd_exclusion_filter;
277 gmx_uint64_t *simd_exclusion_filter64; //!< Used for double w/o SIMD int32 logical support
278 real *simd_interaction_array; /* Array of masks needed for exclusions */
279 int nout; /* The number of force arrays */
280 nbnxn_atomdata_output_t *out; /* Output data structures */
281 int nalloc; /* Allocation size of all arrays (for x/f *x/fstride) */
282 gmx_bool bUseBufferFlags; /* Use the flags or operate on all atoms */
283 nbnxn_buffer_flags_t buffer_flags; /* Flags for buffer zeroing+reduc. */
284 gmx_bool bUseTreeReduce; /* Use tree for force reduction */
285 tMPI_Atomic *syncStep; /* Synchronization step for tree reduce */