e9d47642285cec71163a1a2bc03ce87f127b9efb
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_pairlist.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 #ifndef _nbnxn_pairlist_h
37 #define _nbnxn_pairlist_h
38
39 #include "config.h"
40
41 #include <cstddef>
42
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"
48
49 struct tMPI_Atomic;
50
51 /*! \cond INTERNAL */
52
53 /*! \brief The setup for generating and pruning the nbnxn pair list.
54  *
55  * Without dynamic pruning rlistOuter=rlistInner.
56  */
57 struct NbnxnListParameters
58 {
59     /*! \brief Constructor producing a struct with dynamic pruning disabled
60      */
61     NbnxnListParameters(real rlist) :
62         useDynamicPruning(false),
63         nstlistPrune(-1),
64         rlistOuter(rlist),
65         rlistInner(rlist),
66         numRollingParts(1)
67     {
68     }
69
70     bool useDynamicPruning; //!< Are we using dynamic pair-list pruning
71     int  nstlistPrune;      //!< Pair-list dynamic pruning interval
72     real rlistOuter;        //!< Cut-off of the larger, outer pair-list
73     real rlistInner;        //!< Cut-off of the smaller, inner pair-list
74     int  numRollingParts;   //!< The number parts to divide the pair-list into for rolling pruning, a value of 1 gives no rolling pruning
75 };
76
77 /*! \endcond */
78
79 /* With GPU kernels the i and j cluster size is 8 atoms for CUDA and can be set at compile time for OpenCL */
80 #if GMX_GPU == GMX_GPU_OPENCL
81 static constexpr int c_nbnxnGpuClusterSize = GMX_OCL_NB_CLUSTER_SIZE;
82 #else
83 static constexpr int c_nbnxnGpuClusterSize = 8;
84 #endif
85
86 /* The number of clusters in a super-cluster, used for GPU */
87 static constexpr int c_nbnxnGpuNumClusterPerSupercluster = 8;
88
89 /* With GPU kernels we group cluster pairs in 4 to optimize memory usage
90  * of integers containing 32 bits.
91  */
92 static constexpr int c_nbnxnGpuJgroupSize = 32/c_nbnxnGpuNumClusterPerSupercluster;
93
94 /* In CUDA the number of threads in a warp is 32 and we have cluster pairs
95  * of 8*8=64 atoms, so it's convenient to store data for cluster pair halves.
96  */
97 static constexpr int c_nbnxnGpuClusterpairSplit = 2;
98
99 /* The fixed size of the exclusion mask array for a half cluster pair */
100 static constexpr int c_nbnxnGpuExclSize = c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
101
102 /* A buffer data structure of 64 bytes
103  * to be placed at the beginning and end of structs
104  * to avoid cache invalidation of the real contents
105  * of the struct by writes to neighboring memory.
106  */
107 typedef struct {
108     int dummy[16];
109 } gmx_cache_protect_t;
110
111 /* Abstract type for pair searching data */
112 typedef struct nbnxn_search * nbnxn_search_t;
113
114 /* Function that should return a pointer *ptr to memory
115  * of size nbytes.
116  * Error handling should be done within this function.
117  */
118 typedef void nbnxn_alloc_t (void **ptr, size_t nbytes);
119
120 /* Function that should free the memory pointed to by *ptr.
121  * NULL should not be passed to this function.
122  */
123 typedef void nbnxn_free_t (void *ptr);
124
125 /* This is the actual cluster-pair list j-entry.
126  * cj is the j-cluster.
127  * The interaction bits in excl are indexed i-major, j-minor.
128  * The cj entries are sorted such that ones with exclusions come first.
129  * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
130  * is found, all subsequent j-entries in the i-entry also have full masks.
131  */
132 typedef struct {
133     int          cj;    /* The j-cluster                    */
134     unsigned int excl;  /* The exclusion (interaction) bits */
135 } nbnxn_cj_t;
136
137 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
138  * The upper bits contain information for non-bonded kernel optimization.
139  * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
140  * But three flags can be used to skip interactions, currently only for subc=0
141  * !(shift & NBNXN_CI_DO_LJ(subc))   => we can skip LJ for all pairs
142  * shift & NBNXN_CI_HALF_LJ(subc)    => we can skip LJ for the second half of i
143  * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
144  */
145 #define NBNXN_CI_SHIFT          127
146 #define NBNXN_CI_DO_LJ(subc)    (1<<(7+3*(subc)))
147 #define NBNXN_CI_HALF_LJ(subc)  (1<<(8+3*(subc)))
148 #define NBNXN_CI_DO_COUL(subc)  (1<<(9+3*(subc)))
149
150 /* Simple pair-list i-unit */
151 typedef struct {
152     int ci;             /* i-cluster             */
153     int shift;          /* Shift vector index plus possible flags, see above */
154     int cj_ind_start;   /* Start index into cj   */
155     int cj_ind_end;     /* End index into cj     */
156 } nbnxn_ci_t;
157
158 /* Grouped pair-list i-unit */
159 typedef struct {
160     int sci;            /* i-super-cluster       */
161     int shift;          /* Shift vector index plus possible flags */
162     int cj4_ind_start;  /* Start index into cj4  */
163     int cj4_ind_end;    /* End index into cj4    */
164 } nbnxn_sci_t;
165
166 typedef struct {
167     unsigned int imask;    /* The i-cluster interactions mask for 1 warp  */
168     int          excl_ind; /* Index into the exclusion array for 1 warp   */
169 } nbnxn_im_ei_t;
170
171 typedef struct {
172     int           cj[c_nbnxnGpuJgroupSize];         /* The 4 j-clusters */
173     nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data       for 2 warps   */
174 } nbnxn_cj4_t;
175
176 typedef struct {
177     unsigned int pair[c_nbnxnGpuExclSize]; /* Topology exclusion interaction bits for one warp,
178                                             * each unsigned has bitS for 4*8 i clusters
179                                             */
180 } nbnxn_excl_t;
181
182 typedef struct nbnxn_pairlist_t {
183     gmx_cache_protect_t cp0;
184
185     nbnxn_alloc_t      *alloc;
186     nbnxn_free_t       *free;
187
188     gmx_bool            bSimple;         /* Simple list has na_sc=na_s and uses cj   *
189                                           * Complex list uses cj4                    */
190
191     int                     na_ci;       /* The number of atoms per i-cluster        */
192     int                     na_cj;       /* The number of atoms per j-cluster        */
193     int                     na_sc;       /* The number of atoms per super cluster    */
194     real                    rlist;       /* The radius for constructing the list     */
195     int                     nci;         /* The number of i-clusters in the list     */
196     int                     nciOuter;    /* The number of i-clusters in the outer, unpruned list, -1 when invalid */
197     nbnxn_ci_t             *ci;          /* The i-cluster list, size nci             */
198     nbnxn_ci_t             *ciOuter;     /* The outer, unpruned i-cluster list, size nciOuter(=-1 when invalid) */
199     int                     ci_nalloc;   /* The allocation size of ci/ciOuter        */
200     int                     nsci;        /* The number of i-super-clusters in the list */
201     nbnxn_sci_t            *sci;         /* The i-super-cluster list                 */
202     int                     sci_nalloc;  /* The allocation size of sci               */
203
204     int                     ncj;         /* The number of j-clusters in the list     */
205     nbnxn_cj_t             *cj;          /* The j-cluster list, size ncj             */
206     nbnxn_cj_t             *cjOuter;     /* The outer, unpruned j-cluster list, size ncj    */
207     int                     cj_nalloc;   /* The allocation size of cj/cj0            */
208     int                     ncjInUse;    /* The number of j-clusters that are used by ci entries in this list, will be <= ncj */
209
210     int                     ncj4;        /* The total number of 4*j clusters         */
211     nbnxn_cj4_t            *cj4;         /* The 4*j cluster list, size ncj4          */
212     int                     cj4_nalloc;  /* The allocation size of cj4               */
213     int                     nexcl;       /* The count for excl                       */
214     nbnxn_excl_t           *excl;        /* Atom interaction bits (non-exclusions)   */
215     int                     excl_nalloc; /* The allocation size for excl             */
216     int                     nci_tot;     /* The total number of i clusters           */
217
218     struct nbnxn_list_work *work;
219
220     gmx_cache_protect_t     cp1;
221 } nbnxn_pairlist_t;
222
223 typedef struct {
224     int                nnbl;                  /* number of lists */
225     nbnxn_pairlist_t **nbl;                   /* lists */
226     nbnxn_pairlist_t **nbl_work;              /* work space for rebalancing lists */
227     gmx_bool           bCombined;             /* TRUE if lists get combined into one (the 1st) */
228     gmx_bool           bSimple;               /* TRUE if the list of of type "simple"
229                                                  (na_sc=na_s, no super-clusters used) */
230     int                natpair_ljq;           /* Total number of atom pairs for LJ+Q kernel */
231     int                natpair_lj;            /* Total number of atom pairs for LJ kernel   */
232     int                natpair_q;             /* Total number of atom pairs for Q kernel    */
233     t_nblist         **nbl_fep;               /* List of free-energy atom pair interactions */
234     gmx_int64_t        outerListCreationStep; /* Step at which the outer list was created */
235 } nbnxn_pairlist_set_t;
236
237 enum {
238     nbatXYZ, nbatXYZQ, nbatX4, nbatX8
239 };
240
241 typedef struct {
242     real *f;      /* f, size natoms*fstride                             */
243     real *fshift; /* Shift force array, size SHIFTS*DIM                 */
244     int   nV;     /* The size of *Vvdw and *Vc                          */
245     real *Vvdw;   /* Temporary Van der Waals group energy storage       */
246     real *Vc;     /* Temporary Coulomb group energy storage             */
247     int   nVS;    /* The size of *VSvdw and *VSc                        */
248     real *VSvdw;  /* Temporary SIMD Van der Waals group energy storage  */
249     real *VSc;    /* Temporary SIMD Coulomb group energy storage        */
250 } nbnxn_atomdata_output_t;
251
252 /* Block size in atoms for the non-bonded thread force-buffer reduction,
253  * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
254  * Should be small to reduce the reduction and zeroing cost,
255  * but too small will result in overhead.
256  * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
257  */
258 #if GMX_DOUBLE
259 #define NBNXN_BUFFERFLAG_SIZE   8
260 #else
261 #define NBNXN_BUFFERFLAG_SIZE  16
262 #endif
263
264 /* We store the reduction flags as gmx_bitmask_t.
265  * This limits the number of flags to BITMASK_SIZE.
266  */
267 #define NBNXN_BUFFERFLAG_MAX_THREADS  (BITMASK_SIZE)
268
269 /* Flags for telling if threads write to force output buffers */
270 typedef struct {
271     int               nflag;       /* The number of flag blocks                         */
272     gmx_bitmask_t    *flag;        /* Bit i is set when thread i writes to a cell-block */
273     int               flag_nalloc; /* Allocation size of cxy_flag                       */
274 } nbnxn_buffer_flags_t;
275
276 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
277 enum {
278     ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
279 };
280
281 typedef struct nbnxn_atomdata_t {
282     nbnxn_alloc_t           *alloc;
283     nbnxn_free_t            *free;
284     int                      ntype;           /* The number of different atom types                 */
285     real                    *nbfp;            /* Lennard-Jones 6*C6 and 12*C12 params, size ntype^2*2 */
286     int                      comb_rule;       /* Combination rule, see enum above                   */
287     real                    *nbfp_comb;       /* LJ parameter per atom type, size ntype*2           */
288     real                    *nbfp_aligned;    /* As nbfp, but with an alignment (stride) suitable
289                                                * for the present SIMD architectures
290                                                */
291     int                      natoms;          /* Number of atoms                                    */
292     int                      natoms_local;    /* Number of local atoms                           */
293     int                     *type;            /* Atom types                                         */
294     real                    *lj_comb;         /* LJ parameters per atom for combining for pairs     */
295     int                      XFormat;         /* The format of x (and q), enum                      */
296     int                      FFormat;         /* The format of f, enum                              */
297     real                    *q;               /* Charges, can be NULL if incorporated in x          */
298     int                      na_c;            /* The number of atoms per cluster                    */
299     int                      nenergrp;        /* The number of energy groups                        */
300     int                      neg_2log;        /* Log2 of nenergrp                                   */
301     int                     *energrp;         /* The energy groups per cluster, can be NULL         */
302     gmx_bool                 bDynamicBox;     /* Do we need to update shift_vec every step?    */
303     rvec                    *shift_vec;       /* Shift vectors, copied from t_forcerec              */
304     int                      xstride;         /* stride for a coordinate in x (usually 3 or 4)      */
305     int                      fstride;         /* stride for a coordinate in f (usually 3 or 4)      */
306     real                    *x;               /* x and possibly q, size natoms*xstride              */
307
308     /* j-atom minus i-atom index for generating self and Newton exclusions
309      * cluster-cluster pairs of the diagonal, for 4xn and 2xnn kernels.
310      */
311     real                    *simd_4xn_diagonal_j_minus_i;
312     real                    *simd_2xnn_diagonal_j_minus_i;
313     /* Filters for topology exclusion masks for the SIMD kernels. */
314     gmx_uint32_t            *simd_exclusion_filter;
315     gmx_uint64_t            *simd_exclusion_filter64; //!< Used for double w/o SIMD int32 logical support
316     real                    *simd_interaction_array;  /* Array of masks needed for exclusions */
317     int                      nout;                    /* The number of force arrays                         */
318     nbnxn_atomdata_output_t *out;                     /* Output data structures               */
319     int                      nalloc;                  /* Allocation size of all arrays (for x/f *x/fstride) */
320     gmx_bool                 bUseBufferFlags;         /* Use the flags or operate on all atoms     */
321     nbnxn_buffer_flags_t     buffer_flags;            /* Flags for buffer zeroing+reduc.  */
322     gmx_bool                 bUseTreeReduce;          /* Use tree for force reduction */
323     tMPI_Atomic             *syncStep;                /* Synchronization step for tree reduce */
324 } nbnxn_atomdata_t;
325
326 #endif