Extract nbnxm PairlistSets
[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 NbnxnPairlistCpuWork;
58 struct NbnxnPairlistGpuWork;
59 struct nbnxn_atomdata_t;
60 struct PairlistParams;
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 /* With CPU kernels the i-cluster size is always 4 atoms. */
79 static constexpr int c_nbnxnCpuIClusterSize = 4;
80
81 /* With GPU kernels the i and j cluster size is 8 atoms for CUDA and can be set at compile time for OpenCL */
82 #if GMX_GPU == GMX_GPU_OPENCL
83 static constexpr int c_nbnxnGpuClusterSize = GMX_OPENCL_NB_CLUSTER_SIZE;
84 #else
85 static constexpr int c_nbnxnGpuClusterSize = 8;
86 #endif
87
88 /* The number of clusters in a pair-search cell, used for GPU */
89 static constexpr int c_gpuNumClusterPerCellZ = 2;
90 static constexpr int c_gpuNumClusterPerCellY = 2;
91 static constexpr int c_gpuNumClusterPerCellX = 2;
92 static constexpr int c_gpuNumClusterPerCell  = c_gpuNumClusterPerCellZ*c_gpuNumClusterPerCellY*c_gpuNumClusterPerCellX;
93
94
95 /* In CUDA the number of threads in a warp is 32 and we have cluster pairs
96  * of 8*8=64 atoms, so it's convenient to store data for cluster pair halves.
97  */
98 static constexpr int c_nbnxnGpuClusterpairSplit = 2;
99
100 /* The fixed size of the exclusion mask array for a half cluster pair */
101 static constexpr int c_nbnxnGpuExclSize = c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
102
103 //! The available pair list types
104 enum class PairlistType : int
105 {
106     Simple4x2,
107     Simple4x4,
108     Simple4x8,
109     HierarchicalNxN,
110     Count
111 };
112
113 //! Gives the i-cluster size for each pairlist type
114 static constexpr gmx::EnumerationArray<PairlistType, int> IClusterSizePerListType =
115 { {
116       c_nbnxnCpuIClusterSize,
117       c_nbnxnCpuIClusterSize,
118       c_nbnxnCpuIClusterSize,
119       c_nbnxnGpuClusterSize
120   } };
121 //! Gives the j-cluster size for each pairlist type
122 static constexpr gmx::EnumerationArray<PairlistType, int> JClusterSizePerListType =
123 { {
124       2,
125       4,
126       8,
127       c_nbnxnGpuClusterSize
128   } };
129
130 /* A buffer data structure of 64 bytes
131  * to be placed at the beginning and end of structs
132  * to avoid cache invalidation of the real contents
133  * of the struct by writes to neighboring memory.
134  */
135 typedef struct {
136     int dummy[16];
137 } gmx_cache_protect_t;
138
139 /* This is the actual cluster-pair list j-entry.
140  * cj is the j-cluster.
141  * The interaction bits in excl are indexed i-major, j-minor.
142  * The cj entries are sorted such that ones with exclusions come first.
143  * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
144  * is found, all subsequent j-entries in the i-entry also have full masks.
145  */
146 struct nbnxn_cj_t
147 {
148     int          cj;    /* The j-cluster                    */
149     unsigned int excl;  /* The exclusion (interaction) bits */
150 };
151
152 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
153  * The upper bits contain information for non-bonded kernel optimization.
154  * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
155  * But three flags can be used to skip interactions, currently only for subc=0
156  * !(shift & NBNXN_CI_DO_LJ(subc))   => we can skip LJ for all pairs
157  * shift & NBNXN_CI_HALF_LJ(subc)    => we can skip LJ for the second half of i
158  * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
159  */
160 #define NBNXN_CI_SHIFT          127
161 #define NBNXN_CI_DO_LJ(subc)    (1<<(7+3*(subc)))
162 #define NBNXN_CI_HALF_LJ(subc)  (1<<(8+3*(subc)))
163 #define NBNXN_CI_DO_COUL(subc)  (1<<(9+3*(subc)))
164
165 /* Cluster-pair Interaction masks
166  * Bit i*j-cluster-size + j tells if atom i and j interact.
167  */
168 // TODO: Rename according to convention when moving into Nbnxn namespace
169 /* All interaction mask is the same for all kernels */
170 constexpr unsigned int NBNXN_INTERACTION_MASK_ALL       = 0xffffffffU;
171 /* 4x4 kernel diagonal mask */
172 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG      = 0x08ceU;
173 /* 4x2 kernel diagonal masks */
174 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_0 = 0x0002U;
175 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_1 = 0x002fU;
176 /* 4x8 kernel diagonal masks */
177 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_0 = 0xf0f8fcfeU;
178 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_1 = 0x0080c0e0U;
179
180 /* Simple pair-list i-unit */
181 struct nbnxn_ci_t
182 {
183     int ci;             /* i-cluster             */
184     int shift;          /* Shift vector index plus possible flags, see above */
185     int cj_ind_start;   /* Start index into cj   */
186     int cj_ind_end;     /* End index into cj     */
187 };
188
189 /* Grouped pair-list i-unit */
190 typedef struct {
191     /* Returns the number of j-cluster groups in this entry */
192     int numJClusterGroups() const
193     {
194         return cj4_ind_end - cj4_ind_start;
195     }
196
197     int sci;            /* i-super-cluster       */
198     int shift;          /* Shift vector index plus possible flags */
199     int cj4_ind_start;  /* Start index into cj4  */
200     int cj4_ind_end;    /* End index into cj4    */
201 } nbnxn_sci_t;
202
203 /* Interaction data for a j-group for one warp */
204 struct nbnxn_im_ei_t
205 {
206     // The i-cluster interactions mask for 1 warp
207     unsigned int imask    = 0U;
208     // Index into the exclusion array for 1 warp, default index 0 which means no exclusions
209     int          excl_ind = 0;
210 };
211
212 typedef struct {
213     int           cj[c_nbnxnGpuJgroupSize];         /* The 4 j-clusters */
214     nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data       for 2 warps   */
215 } nbnxn_cj4_t;
216
217 /* Struct for storing the atom-pair interaction bits for a cluster pair in a GPU pairlist */
218 struct nbnxn_excl_t
219 {
220     /* Constructor, sets no exclusions, so all atom pairs interacting */
221     nbnxn_excl_t()
222     {
223         for (unsigned int &pairEntry : pair)
224         {
225             pairEntry = NBNXN_INTERACTION_MASK_ALL;
226         }
227     }
228
229     /* Topology exclusion interaction bits per warp */
230     unsigned int pair[c_nbnxnGpuExclSize];
231 };
232
233 /* Cluster pairlist type for use on CPUs */
234 struct NbnxnPairlistCpu
235 {
236     NbnxnPairlistCpu();
237
238     gmx_cache_protect_t     cp0;
239
240     int                     na_ci;       /* The number of atoms per i-cluster        */
241     int                     na_cj;       /* The number of atoms per j-cluster        */
242     real                    rlist;       /* The radius for constructing the list     */
243     FastVector<nbnxn_ci_t>  ci;          /* The i-cluster list                       */
244     FastVector<nbnxn_ci_t>  ciOuter;     /* The outer, unpruned i-cluster list       */
245
246     FastVector<nbnxn_cj_t>  cj;          /* The j-cluster list, size ncj             */
247     FastVector<nbnxn_cj_t>  cjOuter;     /* The outer, unpruned j-cluster list       */
248     int                     ncjInUse;    /* The number of j-clusters that are used by ci entries in this list, will be <= cj.size() */
249
250     int                     nci_tot;     /* The total number of i clusters           */
251
252     /* Working data storage for list construction */
253     std::unique_ptr<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     /* Working data storage for list construction */
288     std::unique_ptr<NbnxnPairlistGpuWork> work;
289
290     gmx_cache_protect_t                   cp1;
291 };
292
293 /*! \internal
294  * \brief An object that holds the local or non-local pairlists
295  */
296 class PairlistSet
297 {
298     public:
299         //! Constructor: initializes the pairlist set as empty
300         PairlistSet(Nbnxm::InteractionLocality  locality,
301                     const PairlistParams       &listParams);
302
303         ~PairlistSet();
304
305         //! Constructs the pairlists in the set using the coordinates in \p nbat
306         void constructPairlists(PairSearch        *pairSearch,
307                                 nbnxn_atomdata_t  *nbat,
308                                 const t_blocka    *excl,
309                                 Nbnxm::KernelType  kernelType,
310                                 int                minimumIlistCountForGpuBalancing,
311                                 t_nrnb            *nrnb);
312
313         //! Dispatch the kernel for dynamic pairlist pruning
314         void dispatchPruneKernel(const nbnxn_atomdata_t *nbat,
315                                  const rvec             *shift_vec,
316                                  Nbnxm::KernelType       kernelType);
317
318         //! Returns the locality
319         Nbnxm::InteractionLocality locality() const
320         {
321             return locality_;
322         }
323
324         //! Returns the lists of CPU pairlists
325         gmx::ArrayRef<const NbnxnPairlistCpu> cpuLists() const
326         {
327             return cpuLists_;
328         }
329
330         //! Returns a pointer to the GPU pairlist, nullptr when not present
331         const NbnxnPairlistGpu *gpuList() const
332         {
333             if (!gpuLists_.empty())
334             {
335                 return &gpuLists_[0];
336             }
337             else
338             {
339                 return nullptr;
340             }
341         }
342
343         //! Returns the lists of free-energy pairlists, empty when nonbonded interactions are not perturbed
344         gmx::ArrayRef<t_nblist const * const> fepLists() const
345         {
346             return fepLists_;
347         }
348
349     private:
350         //! The locality of the pairlist set
351         Nbnxm::InteractionLocality     locality_;
352         //! List of pairlists in CPU layout
353         std::vector<NbnxnPairlistCpu>  cpuLists_;
354         //! List of working list for rebalancing CPU lists
355         std::vector<NbnxnPairlistCpu>  cpuListsWork_;
356         //! List of pairlists in GPU layout
357         std::vector<NbnxnPairlistGpu>  gpuLists_;
358         //! Pairlist parameters describing setup and ranges
359         const PairlistParams          &params_;
360         //! Tells whether multiple lists get merged into one (the first) after creation
361         bool                           combineLists_;
362         //! Tells whether the lists is of CPU type, otherwise GPU type
363         gmx_bool                       isCpuType_;
364         //! Lists for perturbed interactions in simple atom-atom layout
365         std::vector<t_nblist *>        fepLists_;
366
367     public:
368         /* Pair counts for flop counting */
369         //! Total number of atom pairs for LJ+Q kernel
370         int natpair_ljq_;
371         //! Total number of atom pairs for LJ kernel
372         int natpair_lj_;
373         //! Total number of atom pairs for Q kernel
374         int natpair_q_;
375 };
376
377 //! Initializes a free-energy pair-list
378 void nbnxn_init_pairlist_fep(t_nblist *nl);
379
380 #endif