Convert nbnxn_pairlist_set_t to class PairlistSet
[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/defaultinitializationallocator.h"
48 #include "gromacs/utility/enumerationhelpers.h"
49 #include "gromacs/utility/real.h"
50
51 // This file with constants is separate from this file to be able
52 // to include it during OpenCL jitting without including config.h
53 #include "gromacs/nbnxm/constants.h"
54
55 #include "locality.h"
56
57 struct NbnxnListParameters;
58 struct NbnxnPairlistCpuWork;
59 struct NbnxnPairlistGpuWork;
60 struct nbnxn_atomdata_t;
61 class PairSearch;
62 struct t_blocka;
63 struct t_nrnb;
64
65 namespace Nbnxm
66 {
67 enum class KernelType;
68 }
69
70 /* Convenience type for vector with aligned memory */
71 template<typename T>
72 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
73
74 /* Convenience type for vector that avoids initialization at resize() */
75 template<typename T>
76 using FastVector = std::vector < T, gmx::DefaultInitializationAllocator < T>>;
77
78 enum class PairlistType : int
79 {
80     Simple4x2,
81     Simple4x4,
82     Simple4x8,
83     Hierarchical8x8,
84     Count
85 };
86
87 static constexpr gmx::EnumerationArray<PairlistType, int> IClusterSizePerListType = { { 4, 4, 4, 8 } };
88 static constexpr gmx::EnumerationArray<PairlistType, int> JClusterSizePerListType = { { 2, 4, 8, 8 } };
89
90 /* With CPU kernels the i-cluster size is always 4 atoms. */
91 static constexpr int c_nbnxnCpuIClusterSize = 4;
92
93 /* With GPU kernels the i and j cluster size is 8 atoms for CUDA and can be set at compile time for OpenCL */
94 #if GMX_GPU == GMX_GPU_OPENCL
95 static constexpr int c_nbnxnGpuClusterSize = GMX_OPENCL_NB_CLUSTER_SIZE;
96 #else
97 static constexpr int c_nbnxnGpuClusterSize = 8;
98 #endif
99
100 /* The number of clusters in a pair-search cell, used for GPU */
101 static constexpr int c_gpuNumClusterPerCellZ = 2;
102 static constexpr int c_gpuNumClusterPerCellY = 2;
103 static constexpr int c_gpuNumClusterPerCellX = 2;
104 static constexpr int c_gpuNumClusterPerCell  = c_gpuNumClusterPerCellZ*c_gpuNumClusterPerCellY*c_gpuNumClusterPerCellX;
105
106
107 /* In CUDA the number of threads in a warp is 32 and we have cluster pairs
108  * of 8*8=64 atoms, so it's convenient to store data for cluster pair halves.
109  */
110 static constexpr int c_nbnxnGpuClusterpairSplit = 2;
111
112 /* The fixed size of the exclusion mask array for a half cluster pair */
113 static constexpr int c_nbnxnGpuExclSize = c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
114
115 /* A buffer data structure of 64 bytes
116  * to be placed at the beginning and end of structs
117  * to avoid cache invalidation of the real contents
118  * of the struct by writes to neighboring memory.
119  */
120 typedef struct {
121     int dummy[16];
122 } gmx_cache_protect_t;
123
124 /* This is the actual cluster-pair list j-entry.
125  * cj is the j-cluster.
126  * The interaction bits in excl are indexed i-major, j-minor.
127  * The cj entries are sorted such that ones with exclusions come first.
128  * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
129  * is found, all subsequent j-entries in the i-entry also have full masks.
130  */
131 struct nbnxn_cj_t
132 {
133     int          cj;    /* The j-cluster                    */
134     unsigned int excl;  /* The exclusion (interaction) bits */
135 };
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 /* Cluster-pair Interaction masks
151  * Bit i*j-cluster-size + j tells if atom i and j interact.
152  */
153 // TODO: Rename according to convention when moving into Nbnxn namespace
154 /* All interaction mask is the same for all kernels */
155 constexpr unsigned int NBNXN_INTERACTION_MASK_ALL       = 0xffffffffU;
156 /* 4x4 kernel diagonal mask */
157 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG      = 0x08ceU;
158 /* 4x2 kernel diagonal masks */
159 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_0 = 0x0002U;
160 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_1 = 0x002fU;
161 /* 4x8 kernel diagonal masks */
162 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_0 = 0xf0f8fcfeU;
163 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_1 = 0x0080c0e0U;
164
165 /* Simple pair-list i-unit */
166 struct nbnxn_ci_t
167 {
168     int ci;             /* i-cluster             */
169     int shift;          /* Shift vector index plus possible flags, see above */
170     int cj_ind_start;   /* Start index into cj   */
171     int cj_ind_end;     /* End index into cj     */
172 };
173
174 /* Grouped pair-list i-unit */
175 typedef struct {
176     /* Returns the number of j-cluster groups in this entry */
177     int numJClusterGroups() const
178     {
179         return cj4_ind_end - cj4_ind_start;
180     }
181
182     int sci;            /* i-super-cluster       */
183     int shift;          /* Shift vector index plus possible flags */
184     int cj4_ind_start;  /* Start index into cj4  */
185     int cj4_ind_end;    /* End index into cj4    */
186 } nbnxn_sci_t;
187
188 /* Interaction data for a j-group for one warp */
189 struct nbnxn_im_ei_t
190 {
191     // The i-cluster interactions mask for 1 warp
192     unsigned int imask    = 0U;
193     // Index into the exclusion array for 1 warp, default index 0 which means no exclusions
194     int          excl_ind = 0;
195 };
196
197 typedef struct {
198     int           cj[c_nbnxnGpuJgroupSize];         /* The 4 j-clusters */
199     nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data       for 2 warps   */
200 } nbnxn_cj4_t;
201
202 /* Struct for storing the atom-pair interaction bits for a cluster pair in a GPU pairlist */
203 struct nbnxn_excl_t
204 {
205     /* Constructor, sets no exclusions, so all atom pairs interacting */
206     nbnxn_excl_t()
207     {
208         for (unsigned int &pairEntry : pair)
209         {
210             pairEntry = NBNXN_INTERACTION_MASK_ALL;
211         }
212     }
213
214     /* Topology exclusion interaction bits per warp */
215     unsigned int pair[c_nbnxnGpuExclSize];
216 };
217
218 /* Cluster pairlist type for use on CPUs */
219 struct NbnxnPairlistCpu
220 {
221     NbnxnPairlistCpu();
222
223     gmx_cache_protect_t     cp0;
224
225     int                     na_ci;       /* The number of atoms per i-cluster        */
226     int                     na_cj;       /* The number of atoms per j-cluster        */
227     real                    rlist;       /* The radius for constructing the list     */
228     FastVector<nbnxn_ci_t>  ci;          /* The i-cluster list                       */
229     FastVector<nbnxn_ci_t>  ciOuter;     /* The outer, unpruned i-cluster list       */
230
231     FastVector<nbnxn_cj_t>  cj;          /* The j-cluster list, size ncj             */
232     FastVector<nbnxn_cj_t>  cjOuter;     /* The outer, unpruned j-cluster list       */
233     int                     ncjInUse;    /* The number of j-clusters that are used by ci entries in this list, will be <= cj.size() */
234
235     int                     nci_tot;     /* The total number of i clusters           */
236
237     /* Working data storage for list construction */
238     std::unique_ptr<NbnxnPairlistCpuWork> work;
239
240     gmx_cache_protect_t                   cp1;
241 };
242
243 /* Cluster pairlist type, with extra hierarchies, for on the GPU
244  *
245  * NOTE: for better performance when combining lists over threads,
246  *       all vectors should use default initialization. But when
247  *       changing this, excl should be intialized when adding entries.
248  */
249 struct NbnxnPairlistGpu
250 {
251     /* Constructor
252      *
253      * \param[in] pinningPolicy  Sets the pinning policy for all buffers used on the GPU
254      */
255     NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy);
256
257     gmx_cache_protect_t            cp0;
258
259     int                            na_ci; /* The number of atoms per i-cluster        */
260     int                            na_cj; /* The number of atoms per j-cluster        */
261     int                            na_sc; /* The number of atoms per super cluster    */
262     real                           rlist; /* The radius for constructing the list     */
263     // The i-super-cluster list, indexes into cj4;
264     gmx::HostVector<nbnxn_sci_t>   sci;
265     // The list of 4*j-cluster groups
266     gmx::HostVector<nbnxn_cj4_t>   cj4;
267     // Atom interaction bits (non-exclusions)
268     gmx::HostVector<nbnxn_excl_t>  excl;
269     // The total number of i-clusters
270     int                            nci_tot;
271
272     /* Working data storage for list construction */
273     std::unique_ptr<NbnxnPairlistGpuWork> work;
274
275     gmx_cache_protect_t                   cp1;
276 };
277
278 /*! \internal
279  * \brief An object that holds the local or non-local pairlists
280  */
281 class PairlistSet
282 {
283     public:
284         //! Constructor: initializes the pairlist set as empty
285         PairlistSet(Nbnxm::InteractionLocality  locality,
286                     const NbnxnListParameters  &listParams);
287
288         ~PairlistSet();
289
290         //! Constructs the pairlists in the set using the coordinates in \p nbat
291         void constructPairlists(PairSearch        *pairSearch,
292                                 nbnxn_atomdata_t  *nbat,
293                                 const t_blocka    *excl,
294                                 Nbnxm::KernelType  kernelType,
295                                 int                minimumIlistCountForGpuBalancing,
296                                 t_nrnb            *nrnb);
297
298         //! Dispatch the kernel for dynamic pairlist pruning
299         void dispatchPruneKernel(const nbnxn_atomdata_t *nbat,
300                                  const rvec             *shift_vec,
301                                  Nbnxm::KernelType       kernelType);
302
303         //! Returns the locality
304         Nbnxm::InteractionLocality locality() const
305         {
306             return locality_;
307         }
308
309         //! Returns the lists of CPU pairlists
310         gmx::ArrayRef<const NbnxnPairlistCpu> cpuLists() const
311         {
312             return cpuLists_;
313         }
314
315         //! Returns a pointer to the GPU pairlist, nullptr when not present
316         const NbnxnPairlistGpu *gpuList() const
317         {
318             if (!gpuLists_.empty())
319             {
320                 return &gpuLists_[0];
321             }
322             else
323             {
324                 return nullptr;
325             }
326         }
327
328         //! Returns the lists of free-energy pairlists, empty when nonbonded interactions are not perturbed
329         gmx::ArrayRef<t_nblist const * const> fepLists() const
330         {
331             return fepLists_;
332         }
333
334     private:
335         //! The locality of the pairlist set
336         Nbnxm::InteractionLocality     locality_;
337         //! List of pairlists in CPU layout
338         std::vector<NbnxnPairlistCpu>  cpuLists_;
339         //! List of working list for rebalancing CPU lists
340         std::vector<NbnxnPairlistCpu>  cpuListsWork_;
341         //! List of pairlists in GPU layout
342         std::vector<NbnxnPairlistGpu>  gpuLists_;
343         //! Pairlist parameters describing setup and ranges
344         const NbnxnListParameters     &params_;
345         //! Tells whether multiple lists get merged into one (the first) after creation
346         bool                           combineLists_;
347         //! Tells whether the lists is of CPU type, otherwise GPU type
348         gmx_bool                       isCpuType_;
349         //! Lists for perturbed interactions in simple atom-atom layout
350         std::vector<t_nblist *>        fepLists_;
351
352     public:
353         /* Pair counts for flop counting */
354         //! Total number of atom pairs for LJ+Q kernel
355         int natpair_ljq_;
356         //! Total number of atom pairs for LJ kernel
357         int natpair_lj_;
358         //! Total number of atom pairs for Q kernel
359         int natpair_q_;
360 };
361
362 //! Initializes a free-energy pair-list
363 void nbnxn_init_pairlist_fep(t_nblist *nl);
364
365 #endif