Next-generation SIMD module, reference implementation
[gromacs.git] / src / gromacs / mdlib / nbnxn_pairlist.h
blobe629519ca91ed1565c0938d0d8e861eea54882e9
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, 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
39 #include <cstddef>
41 #include "thread_mpi/atomic.h"
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/mdtypes/nblist.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/bitmask.h"
47 #include "gromacs/utility/real.h"
49 /* A buffer data structure of 64 bytes
50 * to be placed at the beginning and end of structs
51 * to avoid cache invalidation of the real contents
52 * of the struct by writes to neighboring memory.
54 typedef struct {
55 int dummy[16];
56 } gmx_cache_protect_t;
58 /* Abstract type for pair searching data */
59 typedef struct nbnxn_search * nbnxn_search_t;
61 /* Function that should return a pointer *ptr to memory
62 * of size nbytes.
63 * Error handling should be done within this function.
65 typedef void nbnxn_alloc_t (void **ptr, size_t nbytes);
67 /* Function that should free the memory pointed to by *ptr.
68 * NULL should not be passed to this function.
70 typedef void nbnxn_free_t (void *ptr);
72 /* This is the actual cluster-pair list j-entry.
73 * cj is the j-cluster.
74 * The interaction bits in excl are indexed i-major, j-minor.
75 * The cj entries are sorted such that ones with exclusions come first.
76 * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
77 * is found, all subsequent j-entries in the i-entry also have full masks.
79 typedef struct {
80 int cj; /* The j-cluster */
81 unsigned int excl; /* The exclusion (interaction) bits */
82 } nbnxn_cj_t;
84 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
85 * The upper bits contain information for non-bonded kernel optimization.
86 * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
87 * But three flags can be used to skip interactions, currently only for subc=0
88 * !(shift & NBNXN_CI_DO_LJ(subc)) => we can skip LJ for all pairs
89 * shift & NBNXN_CI_HALF_LJ(subc) => we can skip LJ for the second half of i
90 * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
92 #define NBNXN_CI_SHIFT 127
93 #define NBNXN_CI_DO_LJ(subc) (1<<(7+3*(subc)))
94 #define NBNXN_CI_HALF_LJ(subc) (1<<(8+3*(subc)))
95 #define NBNXN_CI_DO_COUL(subc) (1<<(9+3*(subc)))
97 /* Simple pair-list i-unit */
98 typedef struct {
99 int ci; /* i-cluster */
100 int shift; /* Shift vector index plus possible flags, see above */
101 int cj_ind_start; /* Start index into cj */
102 int cj_ind_end; /* End index into cj */
103 } nbnxn_ci_t;
105 /* Grouped pair-list i-unit */
106 typedef struct {
107 int sci; /* i-super-cluster */
108 int shift; /* Shift vector index plus possible flags */
109 int cj4_ind_start; /* Start index into cj4 */
110 int cj4_ind_end; /* End index into cj4 */
111 } nbnxn_sci_t;
113 typedef struct {
114 unsigned int imask; /* The i-cluster interactions mask for 1 warp */
115 int excl_ind; /* Index into the exclusion array for 1 warp */
116 } nbnxn_im_ei_t;
118 typedef struct {
119 int cj[4]; /* The 4 j-clusters */
120 nbnxn_im_ei_t imei[2]; /* The i-cluster mask data for 2 warps */
121 } nbnxn_cj4_t;
123 typedef struct {
124 unsigned int pair[32]; /* Topology exclusion interaction bits for one warp,
125 * each unsigned has bitS for 4*8 i clusters
127 } nbnxn_excl_t;
129 typedef struct nbnxn_pairlist_t {
130 gmx_cache_protect_t cp0;
132 nbnxn_alloc_t *alloc;
133 nbnxn_free_t *free;
135 gmx_bool bSimple; /* Simple list has na_sc=na_s and uses cj *
136 * Complex list uses cj4 */
138 int na_ci; /* The number of atoms per i-cluster */
139 int na_cj; /* The number of atoms per j-cluster */
140 int na_sc; /* The number of atoms per super cluster */
141 real rlist; /* The radius for constructing the list */
142 int nci; /* The number of i-clusters in the list */
143 nbnxn_ci_t *ci; /* The i-cluster list, size nci */
144 int ci_nalloc; /* The allocation size of ci */
145 int nsci; /* The number of i-super-clusters in the list */
146 nbnxn_sci_t *sci; /* The i-super-cluster list */
147 int sci_nalloc; /* The allocation size of sci */
149 int ncj; /* The number of j-clusters in the list */
150 nbnxn_cj_t *cj; /* The j-cluster list, size ncj */
151 int cj_nalloc; /* The allocation size of cj */
153 int ncj4; /* The total number of 4*j clusters */
154 nbnxn_cj4_t *cj4; /* The 4*j cluster list, size ncj4 */
155 int cj4_nalloc; /* The allocation size of cj4 */
156 int nexcl; /* The count for excl */
157 nbnxn_excl_t *excl; /* Atom interaction bits (non-exclusions) */
158 int excl_nalloc; /* The allocation size for excl */
159 int nci_tot; /* The total number of i clusters */
161 struct nbnxn_list_work *work;
163 gmx_cache_protect_t cp1;
164 } nbnxn_pairlist_t;
166 typedef struct {
167 int nnbl; /* number of lists */
168 nbnxn_pairlist_t **nbl; /* lists */
169 gmx_bool bCombined; /* TRUE if lists get combined into one (the 1st) */
170 gmx_bool bSimple; /* TRUE if the list of of type "simple"
171 (na_sc=na_s, no super-clusters used) */
172 int natpair_ljq; /* Total number of atom pairs for LJ+Q kernel */
173 int natpair_lj; /* Total number of atom pairs for LJ kernel */
174 int natpair_q; /* Total number of atom pairs for Q kernel */
175 t_nblist **nbl_fep;
176 } nbnxn_pairlist_set_t;
178 enum {
179 nbatXYZ, nbatXYZQ, nbatX4, nbatX8
182 typedef struct {
183 real *f; /* f, size natoms*fstride */
184 real *fshift; /* Shift force array, size SHIFTS*DIM */
185 int nV; /* The size of *Vvdw and *Vc */
186 real *Vvdw; /* Temporary Van der Waals group energy storage */
187 real *Vc; /* Temporary Coulomb group energy storage */
188 int nVS; /* The size of *VSvdw and *VSc */
189 real *VSvdw; /* Temporary SIMD Van der Waals group energy storage */
190 real *VSc; /* Temporary SIMD Coulomb group energy storage */
191 } nbnxn_atomdata_output_t;
193 /* Block size in atoms for the non-bonded thread force-buffer reduction,
194 * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
195 * Should be small to reduce the reduction and zeroing cost,
196 * but too small will result in overhead.
197 * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
199 #ifdef GMX_DOUBLE
200 #define NBNXN_BUFFERFLAG_SIZE 8
201 #else
202 #define NBNXN_BUFFERFLAG_SIZE 16
203 #endif
205 /* We store the reduction flags as gmx_bitmask_t.
206 * This limits the number of flags to BITMASK_SIZE.
208 #define NBNXN_BUFFERFLAG_MAX_THREADS (BITMASK_SIZE)
210 /* Flags for telling if threads write to force output buffers */
211 typedef struct {
212 int nflag; /* The number of flag blocks */
213 gmx_bitmask_t *flag; /* Bit i is set when thread i writes to a cell-block */
214 int flag_nalloc; /* Allocation size of cxy_flag */
215 } nbnxn_buffer_flags_t;
217 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
218 enum {
219 ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
222 typedef struct nbnxn_atomdata_t {
223 nbnxn_alloc_t *alloc;
224 nbnxn_free_t *free;
225 int ntype; /* The number of different atom types */
226 real *nbfp; /* Lennard-Jones 6*C6 and 12*C12 params, size ntype^2*2 */
227 int comb_rule; /* Combination rule, see enum above */
228 real *nbfp_comb; /* LJ parameter per atom type, size ntype*2 */
229 real *nbfp_aligned; /* As nbfp, but with an alignment (stride) suitable
230 * for the present SIMD architectures
232 int natoms; /* Number of atoms */
233 int natoms_local; /* Number of local atoms */
234 int *type; /* Atom types */
235 real *lj_comb; /* LJ parameters per atom for combining for pairs */
236 int XFormat; /* The format of x (and q), enum */
237 int FFormat; /* The format of f, enum */
238 real *q; /* Charges, can be NULL if incorporated in x */
239 int na_c; /* The number of atoms per cluster */
240 int nenergrp; /* The number of energy groups */
241 int neg_2log; /* Log2 of nenergrp */
242 int *energrp; /* The energy groups per cluster, can be NULL */
243 gmx_bool bDynamicBox; /* Do we need to update shift_vec every step? */
244 rvec *shift_vec; /* Shift vectors, copied from t_forcerec */
245 int xstride; /* stride for a coordinate in x (usually 3 or 4) */
246 int fstride; /* stride for a coordinate in f (usually 3 or 4) */
247 real *x; /* x and possibly q, size natoms*xstride */
249 /* j-atom minus i-atom index for generating self and Newton exclusions
250 * cluster-cluster pairs of the diagonal, for 4xn and 2xnn kernels.
252 real *simd_4xn_diagonal_j_minus_i;
253 real *simd_2xnn_diagonal_j_minus_i;
254 /* Filters for topology exclusion masks for the SIMD kernels. */
255 gmx_uint32_t *simd_exclusion_filter;
256 gmx_uint64_t *simd_exclusion_filter64; //!< Used for double w/o SIMD int32 logical support
257 real *simd_interaction_array; /* Array of masks needed for exclusions */
258 int nout; /* The number of force arrays */
259 nbnxn_atomdata_output_t *out; /* Output data structures */
260 int nalloc; /* Allocation size of all arrays (for x/f *x/fstride) */
261 gmx_bool bUseBufferFlags; /* Use the flags or operate on all atoms */
262 nbnxn_buffer_flags_t buffer_flags; /* Flags for buffer zeroing+reduc. */
263 gmx_bool bUseTreeReduce; /* Use tree for force reduction */
264 tMPI_Atomic_t *syncStep; /* Synchronization step for tree reduce */
265 } nbnxn_atomdata_t;
267 #endif