Convert nbnxn_atomdata_t to C++
[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,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 _nbnxn_pairlist_h
37 #define _nbnxn_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/mdlib/nbnxn_consts.h"
46 #include "gromacs/mdtypes/nblist.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/bitmask.h"
49 #include "gromacs/utility/defaultinitializationallocator.h"
50 #include "gromacs/utility/real.h"
51
52 struct NbnxnPairlistCpuWork;
53 struct NbnxnPairlistGpuWork;
54 struct tMPI_Atomic;
55
56 /* Convenience type for vector with aligned memory */
57 template<typename T>
58 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
59
60 /* Convenience type for vector that avoids initialization at resize() */
61 template<typename T>
62 using FastVector = std::vector < T, gmx::DefaultInitializationAllocator < T>>;
63
64 /*! \cond INTERNAL */
65
66 /*! \brief The setup for generating and pruning the nbnxn pair list.
67  *
68  * Without dynamic pruning rlistOuter=rlistInner.
69  */
70 struct NbnxnListParameters
71 {
72     /*! \brief Constructor producing a struct with dynamic pruning disabled
73      */
74     NbnxnListParameters(real rlist) :
75         useDynamicPruning(false),
76         nstlistPrune(-1),
77         rlistOuter(rlist),
78         rlistInner(rlist),
79         numRollingParts(1)
80     {
81     }
82
83     bool useDynamicPruning; //!< Are we using dynamic pair-list pruning
84     int  nstlistPrune;      //!< Pair-list dynamic pruning interval
85     real rlistOuter;        //!< Cut-off of the larger, outer pair-list
86     real rlistInner;        //!< Cut-off of the smaller, inner pair-list
87     int  numRollingParts;   //!< The number parts to divide the pair-list into for rolling pruning, a value of 1 gives no rolling pruning
88 };
89
90 /*! \endcond */
91
92 /* With CPU kernels the i-cluster size is always 4 atoms. */
93 static constexpr int c_nbnxnCpuIClusterSize = 4;
94
95 /* With GPU kernels the i and j cluster size is 8 atoms for CUDA and can be set at compile time for OpenCL */
96 #if GMX_GPU == GMX_GPU_OPENCL
97 static constexpr int c_nbnxnGpuClusterSize = GMX_OPENCL_NB_CLUSTER_SIZE;
98 #else
99 static constexpr int c_nbnxnGpuClusterSize = 8;
100 #endif
101
102 /* The number of clusters in a pair-search cell, used for GPU */
103 static constexpr int c_gpuNumClusterPerCellZ = 2;
104 static constexpr int c_gpuNumClusterPerCellY = 2;
105 static constexpr int c_gpuNumClusterPerCellX = 2;
106 static constexpr int c_gpuNumClusterPerCell  = c_gpuNumClusterPerCellZ*c_gpuNumClusterPerCellY*c_gpuNumClusterPerCellX;
107
108 /* In CUDA the number of threads in a warp is 32 and we have cluster pairs
109  * of 8*8=64 atoms, so it's convenient to store data for cluster pair halves.
110  */
111 static constexpr int c_nbnxnGpuClusterpairSplit = 2;
112
113 /* The fixed size of the exclusion mask array for a half cluster pair */
114 static constexpr int c_nbnxnGpuExclSize = c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
115
116 /* A buffer data structure of 64 bytes
117  * to be placed at the beginning and end of structs
118  * to avoid cache invalidation of the real contents
119  * of the struct by writes to neighboring memory.
120  */
121 typedef struct {
122     int dummy[16];
123 } gmx_cache_protect_t;
124
125 /* Abstract type for pair searching data */
126 typedef struct nbnxn_search * nbnxn_search_t;
127
128 /* Function that should return a pointer *ptr to memory
129  * of size nbytes.
130  * Error handling should be done within this function.
131  */
132 typedef void nbnxn_alloc_t (void **ptr, size_t nbytes);
133
134 /* Function that should free the memory pointed to by *ptr.
135  * NULL should not be passed to this function.
136  */
137 typedef void nbnxn_free_t (void *ptr);
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 /* 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     int sci;            /* i-super-cluster       */
177     int shift;          /* Shift vector index plus possible flags */
178     int cj4_ind_start;  /* Start index into cj4  */
179     int cj4_ind_end;    /* End index into cj4    */
180 } nbnxn_sci_t;
181
182 typedef struct {
183     unsigned int imask;    /* The i-cluster interactions mask for 1 warp  */
184     int          excl_ind; /* Index into the exclusion array for 1 warp   */
185 } nbnxn_im_ei_t;
186
187 typedef struct {
188     int           cj[c_nbnxnGpuJgroupSize];         /* The 4 j-clusters */
189     nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data       for 2 warps   */
190 } nbnxn_cj4_t;
191
192 typedef struct {
193     unsigned int pair[c_nbnxnGpuExclSize]; /* Topology exclusion interaction bits for one warp,
194                                             * each unsigned has bitS for 4*8 i clusters
195                                             */
196 } nbnxn_excl_t;
197
198 struct NbnxnPairlistCpu
199 {
200     gmx_cache_protect_t     cp0;
201
202     int                     na_ci;       /* The number of atoms per i-cluster        */
203     int                     na_cj;       /* The number of atoms per j-cluster        */
204     real                    rlist;       /* The radius for constructing the list     */
205     FastVector<nbnxn_ci_t>  ci;          /* The i-cluster list                       */
206     FastVector<nbnxn_ci_t>  ciOuter;     /* The outer, unpruned i-cluster list       */
207
208     FastVector<nbnxn_cj_t>  cj;          /* The j-cluster list, size ncj             */
209     FastVector<nbnxn_cj_t>  cjOuter;     /* The outer, unpruned j-cluster list       */
210     int                     ncjInUse;    /* The number of j-clusters that are used by ci entries in this list, will be <= cj.size() */
211
212     int                     nci_tot;     /* The total number of i clusters           */
213
214     NbnxnPairlistCpuWork   *work;
215
216     gmx_cache_protect_t     cp1;
217 };
218
219 struct NbnxnPairlistGpu
220 {
221     gmx_cache_protect_t     cp0;
222
223     nbnxn_alloc_t          *alloc;
224     nbnxn_free_t           *free;
225
226     int                     na_ci;       /* The number of atoms per i-cluster        */
227     int                     na_cj;       /* The number of atoms per j-cluster        */
228     int                     na_sc;       /* The number of atoms per super cluster    */
229     real                    rlist;       /* The radius for constructing the list     */
230     int                     nsci;        /* The number of i-super-clusters in the list */
231     nbnxn_sci_t            *sci;         /* The i-super-cluster list                 */
232     int                     sci_nalloc;  /* The allocation size of sci               */
233
234     int                     ncj4;        /* The total number of 4*j clusters         */
235     nbnxn_cj4_t            *cj4;         /* The 4*j cluster list, size ncj4          */
236     int                     cj4_nalloc;  /* The allocation size of cj4               */
237     int                     nexcl;       /* The count for excl                       */
238     nbnxn_excl_t           *excl;        /* Atom interaction bits (non-exclusions)   */
239     int                     excl_nalloc; /* The allocation size for excl             */
240     int                     nci_tot;     /* The total number of i clusters           */
241
242     NbnxnPairlistGpuWork   *work;
243
244     gmx_cache_protect_t     cp1;
245 };
246
247 typedef struct {
248     int                nnbl;                  /* number of lists */
249     NbnxnPairlistCpu **nbl;                   /* lists for CPU */
250     NbnxnPairlistCpu **nbl_work;              /* work space for rebalancing lists */
251     NbnxnPairlistGpu **nblGpu;                /* lists for GPU */
252     gmx_bool           bCombined;             /* TRUE if lists get combined into one (the 1st) */
253     gmx_bool           bSimple;               /* TRUE if the list of of type "simple"
254                                                  (na_sc=na_s, no super-clusters used) */
255     int                natpair_ljq;           /* Total number of atom pairs for LJ+Q kernel */
256     int                natpair_lj;            /* Total number of atom pairs for LJ kernel   */
257     int                natpair_q;             /* Total number of atom pairs for Q kernel    */
258     t_nblist         **nbl_fep;               /* List of free-energy atom pair interactions */
259     int64_t            outerListCreationStep; /* Step at which the outer list was created */
260 } nbnxn_pairlist_set_t;
261
262 enum {
263     nbatXYZ, nbatXYZQ, nbatX4, nbatX8
264 };
265
266 // Struct that holds force and energy output buffers
267 struct nbnxn_atomdata_output_t
268 {
269     /* Constructor
270      *
271      * \param[in] nb_kernel_type          Type of non-bonded kernel
272      * \param[in] numEnergyGroups         The number of energy groups
273      * \param[in] simdEnergyBufferStride  Stride for entries in the energy buffers for SIMD kernels
274      * \param[in] pinningPolicy           Sets the pinning policy for all buffers used on the GPU
275      */
276     nbnxn_atomdata_output_t(int                nb_kernel_type,
277                             int                numEnergyGroups,
278                             int                simdEnergyBUfferStride,
279                             gmx::PinningPolicy pinningPolicy);
280
281     gmx::HostVector<real> f;      // f, size natoms*fstride
282     gmx::HostVector<real> fshift; // Shift force array, size SHIFTS*DIM
283     gmx::HostVector<real> Vvdw;   // Temporary Van der Waals group energy storage
284     gmx::HostVector<real> Vc;     // Temporary Coulomb group energy storage
285     AlignedVector<real>   VSvdw;  // Temporary SIMD Van der Waals group energy storage
286     AlignedVector<real>   VSc;    // Temporary SIMD Coulomb group energy storage
287 };
288
289 /* Block size in atoms for the non-bonded thread force-buffer reduction,
290  * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
291  * Should be small to reduce the reduction and zeroing cost,
292  * but too small will result in overhead.
293  * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
294  */
295 #if GMX_DOUBLE
296 #define NBNXN_BUFFERFLAG_SIZE   8
297 #else
298 #define NBNXN_BUFFERFLAG_SIZE  16
299 #endif
300
301 /* We store the reduction flags as gmx_bitmask_t.
302  * This limits the number of flags to BITMASK_SIZE.
303  */
304 #define NBNXN_BUFFERFLAG_MAX_THREADS  (BITMASK_SIZE)
305
306 /* Flags for telling if threads write to force output buffers */
307 typedef struct {
308     int               nflag;       /* The number of flag blocks                         */
309     gmx_bitmask_t    *flag;        /* Bit i is set when thread i writes to a cell-block */
310     int               flag_nalloc; /* Allocation size of cxy_flag                       */
311 } nbnxn_buffer_flags_t;
312
313 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
314 enum {
315     ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
316 };
317
318 /* Struct that stores atom related data for the nbnxn module
319  *
320  * Note: performance would improve slightly when all std::vector containers
321  *       in this struct would not initialize during resize().
322  */
323 struct nbnxn_atomdata_t
324 {   //NOLINT(clang-analyzer-optin.performance.Padding)
325     struct Params
326     {
327         /* Constructor
328          *
329          * \param[in] pinningPolicy  Sets the pinning policy for all data that might be transfered to a GPU
330          */
331         Params(gmx::PinningPolicy pinningPolicy);
332
333         // The number of different atom types
334         int                   numTypes;
335         // Lennard-Jone 6*C6 and 12*C12 parameters, size numTypes*2*2
336         gmx::HostVector<real> nbfp;
337         // Combination rule, see enum defined above
338         int                   comb_rule;
339         // LJ parameters per atom type, size numTypes*2
340         gmx::HostVector<real> nbfp_comb;
341         // As nbfp, but with a stride for the present SIMD architecture
342         AlignedVector<real>   nbfp_aligned;
343         // Atom types per atom
344         gmx::HostVector<int>  type;
345         // LJ parameters per atom for fast SIMD loading
346         gmx::HostVector<real> lj_comb;
347         // Charges per atom, not set with format nbatXYZQ
348         gmx::HostVector<real> q;
349         // The number of energy groups
350         int                   nenergrp;
351         // 2log(nenergrp)
352         int                   neg_2log;
353         // The energy groups, one int entry per cluster, only set when needed
354         gmx::HostVector<int>  energrp;
355     };
356
357     // Diagonal and topology exclusion helper data for all SIMD kernels
358     struct SimdMasks
359     {
360         SimdMasks();
361
362         // Helper data for setting up diagonal exclusion masks in the SIMD 4xN kernels
363         AlignedVector<real>     diagonal_4xn_j_minus_i;
364         // Helper data for setting up diaginal exclusion masks in the SIMD 2xNN kernels
365         AlignedVector<real>     diagonal_2xnn_j_minus_i;
366         // Filters for topology exclusion masks for the SIMD kernels
367         AlignedVector<uint32_t> exclusion_filter;
368         // Filters for topology exclusion masks for double SIMD kernels without SIMD int32 logical support
369         AlignedVector<uint64_t> exclusion_filter64;
370         // Array of masks needed for exclusions
371         AlignedVector<real>     interaction_array;
372     };
373
374     /* Constructor
375      *
376      * \param[in] pinningPolicy  Sets the pinning policy for all data that might be transfered to a GPU
377      */
378     nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy);
379
380     /* Returns a const reference to the parameters */
381     const Params &params() const
382     {
383         return params_;
384     }
385
386     /* Returns a non-const reference to the parameters */
387     Params &paramsDeprecated()
388     {
389         return params_;
390     }
391
392     /* Returns the current total number of atoms stored */
393     int numAtoms() const
394     {
395         return numAtoms_;
396     }
397
398     /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
399     gmx::ArrayRef<const real> x() const
400     {
401         return x_;
402     }
403
404     /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
405     gmx::ArrayRef<real> x()
406     {
407         return x_;
408     }
409
410     /* Resizes the coordinate buffer and sets the number of atoms */
411     void resizeCoordinateBuffer(int numAtoms);
412
413     /* Resizes the force buffers for the current number of atoms */
414     void resizeForceBuffers();
415
416     private:
417         // The LJ and charge parameters
418         Params                     params_;
419         // The total number of atoms currently stored
420         int                        numAtoms_;
421     public:
422         int                        natoms_local; /* Number of local atoms                           */
423         int                        XFormat;      /* The format of x (and q), enum                      */
424         int                        FFormat;      /* The format of f, enum                              */
425         gmx_bool                   bDynamicBox;  /* Do we need to update shift_vec every step?    */
426         gmx::HostVector<gmx::RVec> shift_vec;    /* Shift vectors, copied from t_forcerec              */
427         int                        xstride;      /* stride for a coordinate in x (usually 3 or 4)      */
428         int                        fstride;      /* stride for a coordinate in f (usually 3 or 4)      */
429     private:
430         gmx::HostVector<real>      x_;           /* x and possibly q, size natoms*xstride              */
431
432     public:
433         // Masks for handling exclusions in the SIMD kernels
434         const SimdMasks          simdMasks;
435
436         /* Output data */
437         std::vector<nbnxn_atomdata_output_t> out; /* Output data structures, 1 per thread */
438
439         /* Reduction related data */
440         gmx_bool                 bUseBufferFlags;     /* Use the flags or operate on all atoms     */
441         nbnxn_buffer_flags_t     buffer_flags;        /* Flags for buffer zeroing+reduc.  */
442         gmx_bool                 bUseTreeReduce;      /* Use tree for force reduction */
443         tMPI_Atomic             *syncStep;            /* Synchronization step for tree reduce */
444 };
445
446 #endif