Minor nbnxm cleanup
[alexxy/gromacs.git] / src / gromacs / nbnxm / 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,2019, 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 GMX_NBNXM_PAIRLIST_H
37 #define GMX_NBNXM_PAIRLIST_H
38
39 #include "config.h"
40
41 #include <cstddef>
42
43 #include "gromacs/gpu_utils/hostallocator.h"
44 #include "gromacs/math/vectypes.h"
45 #include "gromacs/mdtypes/nblist.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/bitmask.h"
48 #include "gromacs/utility/defaultinitializationallocator.h"
49 #include "gromacs/utility/enumerationhelpers.h"
50 #include "gromacs/utility/real.h"
51
52 // This file with constants is separate from this file to be able
53 // to include it during OpenCL jitting without including config.h
54 #include "gromacs/nbnxm/constants.h"
55
56 struct NbnxnPairlistCpuWork;
57 struct NbnxnPairlistGpuWork;
58 struct tMPI_Atomic;
59
60 namespace Nbnxm
61 {
62 enum class KernelType;
63 }
64
65 /* Convenience type for vector with aligned memory */
66 template<typename T>
67 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
68
69 /* Convenience type for vector that avoids initialization at resize() */
70 template<typename T>
71 using FastVector = std::vector < T, gmx::DefaultInitializationAllocator < T>>;
72
73 enum class PairlistType : int
74 {
75     Simple4x2,
76     Simple4x4,
77     Simple4x8,
78     Hierarchical8x8,
79     Count
80 };
81
82 static constexpr gmx::EnumerationArray<PairlistType, int> IClusterSizePerListType = { 4, 4, 4, 8 };
83 static constexpr gmx::EnumerationArray<PairlistType, int> JClusterSizePerListType = { 2, 4, 8, 8 };
84
85 /*! \cond INTERNAL */
86
87 /*! \brief The setup for generating and pruning the nbnxn pair list.
88  *
89  * Without dynamic pruning rlistOuter=rlistInner.
90  */
91 struct NbnxnListParameters
92 {
93     /*! \brief Constructor producing a struct with dynamic pruning disabled
94      */
95     NbnxnListParameters(Nbnxm::KernelType kernelType,
96                         real              rlist);
97
98     PairlistType pairlistType;      //!< The type of cluster-pair list
99     bool         useDynamicPruning; //!< Are we using dynamic pair-list pruning
100     int          nstlistPrune;      //!< Pair-list dynamic pruning interval
101     real         rlistOuter;        //!< Cut-off of the larger, outer pair-list
102     real         rlistInner;        //!< Cut-off of the smaller, inner pair-list
103     int          numRollingParts;   //!< The number parts to divide the pair-list into for rolling pruning, a value of 1 gives no rolling pruning
104 };
105
106 /*! \endcond */
107
108 /* With CPU kernels the i-cluster size is always 4 atoms. */
109 static constexpr int c_nbnxnCpuIClusterSize = 4;
110
111 /* With GPU kernels the i and j cluster size is 8 atoms for CUDA and can be set at compile time for OpenCL */
112 #if GMX_GPU == GMX_GPU_OPENCL
113 static constexpr int c_nbnxnGpuClusterSize = GMX_OPENCL_NB_CLUSTER_SIZE;
114 #else
115 static constexpr int c_nbnxnGpuClusterSize = 8;
116 #endif
117
118 /* The number of clusters in a pair-search cell, used for GPU */
119 static constexpr int c_gpuNumClusterPerCellZ = 2;
120 static constexpr int c_gpuNumClusterPerCellY = 2;
121 static constexpr int c_gpuNumClusterPerCellX = 2;
122 static constexpr int c_gpuNumClusterPerCell  = c_gpuNumClusterPerCellZ*c_gpuNumClusterPerCellY*c_gpuNumClusterPerCellX;
123
124
125 /* In CUDA the number of threads in a warp is 32 and we have cluster pairs
126  * of 8*8=64 atoms, so it's convenient to store data for cluster pair halves.
127  */
128 static constexpr int c_nbnxnGpuClusterpairSplit = 2;
129
130 /* The fixed size of the exclusion mask array for a half cluster pair */
131 static constexpr int c_nbnxnGpuExclSize = c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
132
133 /* A buffer data structure of 64 bytes
134  * to be placed at the beginning and end of structs
135  * to avoid cache invalidation of the real contents
136  * of the struct by writes to neighboring memory.
137  */
138 typedef struct {
139     int dummy[16];
140 } gmx_cache_protect_t;
141
142 /* This is the actual cluster-pair list j-entry.
143  * cj is the j-cluster.
144  * The interaction bits in excl are indexed i-major, j-minor.
145  * The cj entries are sorted such that ones with exclusions come first.
146  * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
147  * is found, all subsequent j-entries in the i-entry also have full masks.
148  */
149 struct nbnxn_cj_t
150 {
151     int          cj;    /* The j-cluster                    */
152     unsigned int excl;  /* The exclusion (interaction) bits */
153 };
154
155 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
156  * The upper bits contain information for non-bonded kernel optimization.
157  * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
158  * But three flags can be used to skip interactions, currently only for subc=0
159  * !(shift & NBNXN_CI_DO_LJ(subc))   => we can skip LJ for all pairs
160  * shift & NBNXN_CI_HALF_LJ(subc)    => we can skip LJ for the second half of i
161  * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
162  */
163 #define NBNXN_CI_SHIFT          127
164 #define NBNXN_CI_DO_LJ(subc)    (1<<(7+3*(subc)))
165 #define NBNXN_CI_HALF_LJ(subc)  (1<<(8+3*(subc)))
166 #define NBNXN_CI_DO_COUL(subc)  (1<<(9+3*(subc)))
167
168 /* Cluster-pair Interaction masks
169  * Bit i*j-cluster-size + j tells if atom i and j interact.
170  */
171 // TODO: Rename according to convention when moving into Nbnxn namespace
172 /* All interaction mask is the same for all kernels */
173 constexpr unsigned int NBNXN_INTERACTION_MASK_ALL       = 0xffffffffU;
174 /* 4x4 kernel diagonal mask */
175 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG      = 0x08ceU;
176 /* 4x2 kernel diagonal masks */
177 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_0 = 0x0002U;
178 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_1 = 0x002fU;
179 /* 4x8 kernel diagonal masks */
180 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_0 = 0xf0f8fcfeU;
181 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_1 = 0x0080c0e0U;
182
183 /* Simple pair-list i-unit */
184 struct nbnxn_ci_t
185 {
186     int ci;             /* i-cluster             */
187     int shift;          /* Shift vector index plus possible flags, see above */
188     int cj_ind_start;   /* Start index into cj   */
189     int cj_ind_end;     /* End index into cj     */
190 };
191
192 /* Grouped pair-list i-unit */
193 typedef struct {
194     /* Returns the number of j-cluster groups in this entry */
195     int numJClusterGroups() const
196     {
197         return cj4_ind_end - cj4_ind_start;
198     }
199
200     int sci;            /* i-super-cluster       */
201     int shift;          /* Shift vector index plus possible flags */
202     int cj4_ind_start;  /* Start index into cj4  */
203     int cj4_ind_end;    /* End index into cj4    */
204 } nbnxn_sci_t;
205
206 /* Interaction data for a j-group for one warp */
207 struct nbnxn_im_ei_t
208 {
209     // The i-cluster interactions mask for 1 warp
210     unsigned int imask    = 0U;
211     // Index into the exclusion array for 1 warp, default index 0 which means no exclusions
212     int          excl_ind = 0;
213 };
214
215 typedef struct {
216     int           cj[c_nbnxnGpuJgroupSize];         /* The 4 j-clusters */
217     nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data       for 2 warps   */
218 } nbnxn_cj4_t;
219
220 /* Struct for storing the atom-pair interaction bits for a cluster pair in a GPU pairlist */
221 struct nbnxn_excl_t
222 {
223     /* Constructor, sets no exclusions, so all atom pairs interacting */
224     nbnxn_excl_t()
225     {
226         for (unsigned int &pairEntry : pair)
227         {
228             pairEntry = NBNXN_INTERACTION_MASK_ALL;
229         }
230     }
231
232     /* Topology exclusion interaction bits per warp */
233     unsigned int pair[c_nbnxnGpuExclSize];
234 };
235
236 /* Cluster pairlist type for use on CPUs */
237 struct NbnxnPairlistCpu
238 {
239     gmx_cache_protect_t     cp0;
240
241     int                     na_ci;       /* The number of atoms per i-cluster        */
242     int                     na_cj;       /* The number of atoms per j-cluster        */
243     real                    rlist;       /* The radius for constructing the list     */
244     FastVector<nbnxn_ci_t>  ci;          /* The i-cluster list                       */
245     FastVector<nbnxn_ci_t>  ciOuter;     /* The outer, unpruned i-cluster list       */
246
247     FastVector<nbnxn_cj_t>  cj;          /* The j-cluster list, size ncj             */
248     FastVector<nbnxn_cj_t>  cjOuter;     /* The outer, unpruned j-cluster list       */
249     int                     ncjInUse;    /* The number of j-clusters that are used by ci entries in this list, will be <= cj.size() */
250
251     int                     nci_tot;     /* The total number of i clusters           */
252
253     NbnxnPairlistCpuWork   *work;
254
255     gmx_cache_protect_t     cp1;
256 };
257
258 /* Cluster pairlist type, with extra hierarchies, for on the GPU
259  *
260  * NOTE: for better performance when combining lists over threads,
261  *       all vectors should use default initialization. But when
262  *       changing this, excl should be intialized when adding entries.
263  */
264 struct NbnxnPairlistGpu
265 {
266     /* Constructor
267      *
268      * \param[in] pinningPolicy  Sets the pinning policy for all buffers used on the GPU
269      */
270     NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy);
271
272     gmx_cache_protect_t            cp0;
273
274     int                            na_ci; /* The number of atoms per i-cluster        */
275     int                            na_cj; /* The number of atoms per j-cluster        */
276     int                            na_sc; /* The number of atoms per super cluster    */
277     real                           rlist; /* The radius for constructing the list     */
278     // The i-super-cluster list, indexes into cj4;
279     gmx::HostVector<nbnxn_sci_t>   sci;
280     // The list of 4*j-cluster groups
281     gmx::HostVector<nbnxn_cj4_t>   cj4;
282     // Atom interaction bits (non-exclusions)
283     gmx::HostVector<nbnxn_excl_t>  excl;
284     // The total number of i-clusters
285     int                            nci_tot;
286
287     NbnxnPairlistGpuWork          *work;
288
289     gmx_cache_protect_t            cp1;
290 };
291
292 struct nbnxn_pairlist_set_t
293 {
294     nbnxn_pairlist_set_t(const NbnxnListParameters &listParams);
295
296     int                         nnbl;         /* number of lists */
297     NbnxnPairlistCpu          **nbl;          /* lists for CPU */
298     NbnxnPairlistCpu          **nbl_work;     /* work space for rebalancing lists */
299     NbnxnPairlistGpu          **nblGpu;       /* lists for GPU */
300     const NbnxnListParameters  &params;       /* Pairlist parameters desribing setup and ranges */
301     gmx_bool                    bCombined;    /* TRUE if lists get combined into one (the 1st) */
302     gmx_bool                    bSimple;      /* TRUE if the list of of type "simple"
303                                                  (na_sc=na_s, no super-clusters used) */
304
305     /* Counts for debug printing */
306     int                     natpair_ljq;           /* Total number of atom pairs for LJ+Q kernel */
307     int                     natpair_lj;            /* Total number of atom pairs for LJ kernel   */
308     int                     natpair_q;             /* Total number of atom pairs for Q kernel    */
309     std::vector<t_nblist *> nbl_fep;               /* List of free-energy atom pair interactions */
310     int64_t                 outerListCreationStep; /* Step at which the outer list was created */
311 };
312
313 enum {
314     nbatXYZ, nbatXYZQ, nbatX4, nbatX8
315 };
316
317 // Struct that holds force and energy output buffers
318 struct nbnxn_atomdata_output_t
319 {
320     /* Constructor
321      *
322      * \param[in] kernelType              Type of non-bonded kernel
323      * \param[in] numEnergyGroups         The number of energy groups
324      * \param[in] simdEnergyBufferStride  Stride for entries in the energy buffers for SIMD kernels
325      * \param[in] pinningPolicy           Sets the pinning policy for all buffers used on the GPU
326      */
327     nbnxn_atomdata_output_t(Nbnxm::KernelType  kernelType,
328                             int                numEnergyGroups,
329                             int                simdEnergyBUfferStride,
330                             gmx::PinningPolicy pinningPolicy);
331
332     gmx::HostVector<real> f;      // f, size natoms*fstride
333     gmx::HostVector<real> fshift; // Shift force array, size SHIFTS*DIM
334     gmx::HostVector<real> Vvdw;   // Temporary Van der Waals group energy storage
335     gmx::HostVector<real> Vc;     // Temporary Coulomb group energy storage
336     AlignedVector<real>   VSvdw;  // Temporary SIMD Van der Waals group energy storage
337     AlignedVector<real>   VSc;    // Temporary SIMD Coulomb group energy storage
338 };
339
340 /* Block size in atoms for the non-bonded thread force-buffer reduction,
341  * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
342  * Should be small to reduce the reduction and zeroing cost,
343  * but too small will result in overhead.
344  * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
345  */
346 #if GMX_DOUBLE
347 #define NBNXN_BUFFERFLAG_SIZE   8
348 #else
349 #define NBNXN_BUFFERFLAG_SIZE  16
350 #endif
351
352 /* We store the reduction flags as gmx_bitmask_t.
353  * This limits the number of flags to BITMASK_SIZE.
354  */
355 #define NBNXN_BUFFERFLAG_MAX_THREADS  (BITMASK_SIZE)
356
357 /* Flags for telling if threads write to force output buffers */
358 typedef struct {
359     int               nflag;       /* The number of flag blocks                         */
360     gmx_bitmask_t    *flag;        /* Bit i is set when thread i writes to a cell-block */
361     int               flag_nalloc; /* Allocation size of cxy_flag                       */
362 } nbnxn_buffer_flags_t;
363
364 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
365 enum {
366     ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
367 };
368
369 /* Struct that stores atom related data for the nbnxn module
370  *
371  * Note: performance would improve slightly when all std::vector containers
372  *       in this struct would not initialize during resize().
373  */
374 struct nbnxn_atomdata_t
375 {   //NOLINT(clang-analyzer-optin.performance.Padding)
376     struct Params
377     {
378         /* Constructor
379          *
380          * \param[in] pinningPolicy  Sets the pinning policy for all data that might be transfered to a GPU
381          */
382         Params(gmx::PinningPolicy pinningPolicy);
383
384         // The number of different atom types
385         int                   numTypes;
386         // Lennard-Jone 6*C6 and 12*C12 parameters, size numTypes*2*2
387         gmx::HostVector<real> nbfp;
388         // Combination rule, see enum defined above
389         int                   comb_rule;
390         // LJ parameters per atom type, size numTypes*2
391         gmx::HostVector<real> nbfp_comb;
392         // As nbfp, but with a stride for the present SIMD architecture
393         AlignedVector<real>   nbfp_aligned;
394         // Atom types per atom
395         gmx::HostVector<int>  type;
396         // LJ parameters per atom for fast SIMD loading
397         gmx::HostVector<real> lj_comb;
398         // Charges per atom, not set with format nbatXYZQ
399         gmx::HostVector<real> q;
400         // The number of energy groups
401         int                   nenergrp;
402         // 2log(nenergrp)
403         int                   neg_2log;
404         // The energy groups, one int entry per cluster, only set when needed
405         gmx::HostVector<int>  energrp;
406     };
407
408     // Diagonal and topology exclusion helper data for all SIMD kernels
409     struct SimdMasks
410     {
411         SimdMasks();
412
413         // Helper data for setting up diagonal exclusion masks in the SIMD 4xN kernels
414         AlignedVector<real>     diagonal_4xn_j_minus_i;
415         // Helper data for setting up diaginal exclusion masks in the SIMD 2xNN kernels
416         AlignedVector<real>     diagonal_2xnn_j_minus_i;
417         // Filters for topology exclusion masks for the SIMD kernels
418         AlignedVector<uint32_t> exclusion_filter;
419         // Filters for topology exclusion masks for double SIMD kernels without SIMD int32 logical support
420         AlignedVector<uint64_t> exclusion_filter64;
421         // Array of masks needed for exclusions
422         AlignedVector<real>     interaction_array;
423     };
424
425     /* Constructor
426      *
427      * \param[in] pinningPolicy  Sets the pinning policy for all data that might be transfered to a GPU
428      */
429     nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy);
430
431     /* Returns a const reference to the parameters */
432     const Params &params() const
433     {
434         return params_;
435     }
436
437     /* Returns a non-const reference to the parameters */
438     Params &paramsDeprecated()
439     {
440         return params_;
441     }
442
443     /* Returns the current total number of atoms stored */
444     int numAtoms() const
445     {
446         return numAtoms_;
447     }
448
449     /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
450     gmx::ArrayRef<const real> x() const
451     {
452         return x_;
453     }
454
455     /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
456     gmx::ArrayRef<real> x()
457     {
458         return x_;
459     }
460
461     /* Resizes the coordinate buffer and sets the number of atoms */
462     void resizeCoordinateBuffer(int numAtoms);
463
464     /* Resizes the force buffers for the current number of atoms */
465     void resizeForceBuffers();
466
467     private:
468         // The LJ and charge parameters
469         Params                     params_;
470         // The total number of atoms currently stored
471         int                        numAtoms_;
472     public:
473         int                        natoms_local; /* Number of local atoms                           */
474         int                        XFormat;      /* The format of x (and q), enum                      */
475         int                        FFormat;      /* The format of f, enum                              */
476         gmx_bool                   bDynamicBox;  /* Do we need to update shift_vec every step?    */
477         gmx::HostVector<gmx::RVec> shift_vec;    /* Shift vectors, copied from t_forcerec              */
478         int                        xstride;      /* stride for a coordinate in x (usually 3 or 4)      */
479         int                        fstride;      /* stride for a coordinate in f (usually 3 or 4)      */
480     private:
481         gmx::HostVector<real>      x_;           /* x and possibly q, size natoms*xstride              */
482
483     public:
484         // Masks for handling exclusions in the SIMD kernels
485         const SimdMasks          simdMasks;
486
487         /* Output data */
488         std::vector<nbnxn_atomdata_output_t> out; /* Output data structures, 1 per thread */
489
490         /* Reduction related data */
491         gmx_bool                 bUseBufferFlags;     /* Use the flags or operate on all atoms     */
492         nbnxn_buffer_flags_t     buffer_flags;        /* Flags for buffer zeroing+reduc.  */
493         gmx_bool                 bUseTreeReduce;      /* Use tree for force reduction */
494         tMPI_Atomic             *syncStep;            /* Synchronization step for tree reduce */
495 };
496
497 #endif