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