Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / nbnxm / atomdata.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_NBNXN_ATOMDATA_H
37 #define GMX_NBNXN_ATOMDATA_H
38
39 #include <cstdio>
40
41 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
42 #include "gromacs/gpu_utils/hostallocator.h"
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/mdtypes/locality.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/bitmask.h"
47 #include "gromacs/utility/real.h"
48
49 #include "gpu_types.h"
50
51 namespace gmx
52 {
53 class MDLogger;
54 }
55
56 struct nbnxn_atomdata_t;
57 struct nonbonded_verlet_t;
58 struct t_mdatoms;
59 struct tMPI_Atomic;
60
61 class GpuEventSynchronizer;
62
63 namespace Nbnxm
64 {
65 class GridSet;
66 enum class KernelType;
67 } // namespace Nbnxm
68
69 /* Convenience type for vector with aligned memory */
70 template<typename T>
71 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
72
73 enum
74 {
75     nbatXYZ,
76     nbatXYZQ,
77     nbatX4,
78     nbatX8
79 };
80
81 //! Stride for coordinate/force arrays with xyz coordinate storage
82 static constexpr int STRIDE_XYZ = 3;
83 //! Stride for coordinate/force arrays with xyzq coordinate storage
84 static constexpr int STRIDE_XYZQ = 4;
85 //! Size of packs of x, y or z with SIMD 4-grouped packed coordinates/forces
86 static constexpr int c_packX4 = 4;
87 //! Size of packs of x, y or z with SIMD 8-grouped packed coordinates/forces
88 static constexpr int c_packX8 = 8;
89 //! Stridefor a pack of 4 coordinates/forces
90 static constexpr int STRIDE_P4 = DIM * c_packX4;
91 //! Stridefor a pack of 8 coordinates/forces
92 static constexpr int STRIDE_P8 = DIM * c_packX8;
93
94 //! Returns the index in a coordinate array corresponding to atom a
95 template<int packSize>
96 static inline int atom_to_x_index(int a)
97 {
98     return DIM * (a & ~(packSize - 1)) + (a & (packSize - 1));
99 }
100
101 // Struct that holds force and energy output buffers
102 struct nbnxn_atomdata_output_t
103 {
104     /* Constructor
105      *
106      * \param[in] kernelType              Type of non-bonded kernel
107      * \param[in] numEnergyGroups         The number of energy groups
108      * \param[in] simdEnergyBufferStride  Stride for entries in the energy buffers for SIMD kernels
109      * \param[in] pinningPolicy           Sets the pinning policy for all buffers used on the GPU
110      */
111     nbnxn_atomdata_output_t(Nbnxm::KernelType  kernelType,
112                             int                numEnergyGroups,
113                             int                simdEnergyBUfferStride,
114                             gmx::PinningPolicy pinningPolicy);
115
116     gmx::HostVector<real> f;      // f, size natoms*fstride
117     gmx::HostVector<real> fshift; // Shift force array, size SHIFTS*DIM
118     gmx::HostVector<real> Vvdw;   // Temporary Van der Waals group energy storage
119     gmx::HostVector<real> Vc;     // Temporary Coulomb group energy storage
120     AlignedVector<real>   VSvdw;  // Temporary SIMD Van der Waals group energy storage
121     AlignedVector<real>   VSc;    // Temporary SIMD Coulomb group energy storage
122 };
123
124 /* Block size in atoms for the non-bonded thread force-buffer reduction,
125  * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
126  * Should be small to reduce the reduction and zeroing cost,
127  * but too small will result in overhead.
128  * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
129  */
130 #if GMX_DOUBLE
131 #    define NBNXN_BUFFERFLAG_SIZE 8
132 #else
133 #    define NBNXN_BUFFERFLAG_SIZE 16
134 #endif
135
136 /* We store the reduction flags as gmx_bitmask_t.
137  * This limits the number of flags to BITMASK_SIZE.
138  */
139 #define NBNXN_BUFFERFLAG_MAX_THREADS (BITMASK_SIZE)
140
141 /* Flags for telling if threads write to force output buffers */
142 typedef struct
143 {
144     int            nflag;       /* The number of flag blocks                         */
145     gmx_bitmask_t* flag;        /* Bit i is set when thread i writes to a cell-block */
146     int            flag_nalloc; /* Allocation size of cxy_flag                       */
147 } nbnxn_buffer_flags_t;
148
149 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
150 enum
151 {
152     ljcrGEOM,
153     ljcrLB,
154     ljcrNONE,
155     ljcrNR
156 };
157
158 /* Struct that stores atom related data for the nbnxn module
159  *
160  * Note: performance would improve slightly when all std::vector containers
161  *       in this struct would not initialize during resize().
162  */
163 struct nbnxn_atomdata_t
164 { //NOLINT(clang-analyzer-optin.performance.Padding)
165     struct Params
166     {
167         /* Constructor
168          *
169          * \param[in] pinningPolicy  Sets the pinning policy for all data that might be transfered to a GPU
170          */
171         Params(gmx::PinningPolicy pinningPolicy);
172
173         // The number of different atom types
174         int numTypes;
175         // Lennard-Jone 6*C6 and 12*C12 parameters, size numTypes*2*2
176         gmx::HostVector<real> nbfp;
177         // Combination rule, see enum defined above
178         int comb_rule;
179         // LJ parameters per atom type, size numTypes*2
180         gmx::HostVector<real> nbfp_comb;
181         // As nbfp, but with a stride for the present SIMD architecture
182         AlignedVector<real> nbfp_aligned;
183         // Atom types per atom
184         gmx::HostVector<int> type;
185         // LJ parameters per atom for fast SIMD loading
186         gmx::HostVector<real> lj_comb;
187         // Charges per atom, not set with format nbatXYZQ
188         gmx::HostVector<real> q;
189         // The number of energy groups
190         int nenergrp;
191         // 2log(nenergrp)
192         int neg_2log;
193         // The energy groups, one int entry per cluster, only set when needed
194         gmx::HostVector<int> energrp;
195     };
196
197     // Diagonal and topology exclusion helper data for all SIMD kernels
198     struct SimdMasks
199     {
200         SimdMasks();
201
202         // Helper data for setting up diagonal exclusion masks in the SIMD 4xN kernels
203         AlignedVector<real> diagonal_4xn_j_minus_i;
204         // Helper data for setting up diaginal exclusion masks in the SIMD 2xNN kernels
205         AlignedVector<real> diagonal_2xnn_j_minus_i;
206         // Filters for topology exclusion masks for the SIMD kernels
207         AlignedVector<uint32_t> exclusion_filter;
208         // Filters for topology exclusion masks for double SIMD kernels without SIMD int32 logical support
209         AlignedVector<uint64_t> exclusion_filter64;
210         // Array of masks needed for exclusions
211         AlignedVector<real> interaction_array;
212     };
213
214     /* Constructor
215      *
216      * \param[in] pinningPolicy  Sets the pinning policy for all data that might be transfered to a GPU
217      */
218     nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy);
219
220     /* Returns a const reference to the parameters */
221     const Params& params() const { return params_; }
222
223     /* Returns a non-const reference to the parameters */
224     Params& paramsDeprecated() { return params_; }
225
226     /* Returns the current total number of atoms stored */
227     int numAtoms() const { return numAtoms_; }
228
229     /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
230     gmx::ArrayRef<const real> x() const { return x_; }
231
232     /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
233     gmx::ArrayRef<real> x() { return x_; }
234
235     /* Resizes the coordinate buffer and sets the number of atoms */
236     void resizeCoordinateBuffer(int numAtoms);
237
238     /* Resizes the force buffers for the current number of atoms */
239     void resizeForceBuffers();
240
241 private:
242     // The LJ and charge parameters
243     Params params_;
244     // The total number of atoms currently stored
245     int numAtoms_;
246
247 public:
248     int                        natoms_local; /* Number of local atoms                           */
249     int                        XFormat;     /* The format of x (and q), enum                      */
250     int                        FFormat;     /* The format of f, enum                              */
251     gmx_bool                   bDynamicBox; /* Do we need to update shift_vec every step?    */
252     gmx::HostVector<gmx::RVec> shift_vec;   /* Shift vectors, copied from t_forcerec              */
253     int                        xstride;     /* stride for a coordinate in x (usually 3 or 4)      */
254     int                        fstride;     /* stride for a coordinate in f (usually 3 or 4)      */
255 private:
256     gmx::HostVector<real> x_; /* x and possibly q, size natoms*xstride              */
257
258 public:
259     // Masks for handling exclusions in the SIMD kernels
260     const SimdMasks simdMasks;
261
262     /* Output data */
263     std::vector<nbnxn_atomdata_output_t> out; /* Output data structures, 1 per thread */
264
265     /* Reduction related data */
266     gmx_bool             bUseBufferFlags; /* Use the flags or operate on all atoms     */
267     nbnxn_buffer_flags_t buffer_flags;    /* Flags for buffer zeroing+reduc.  */
268     gmx_bool             bUseTreeReduce;  /* Use tree for force reduction */
269     tMPI_Atomic*         syncStep;        /* Synchronization step for tree reduce */
270 };
271
272 /* Copy na rvec elements from x to xnb using nbatFormat, start dest a0,
273  * and fills up to na_round with coordinates that are far away.
274  */
275 void copy_rvec_to_nbat_real(const int* a, int na, int na_round, const rvec* x, int nbatFormat, real* xnb, int a0);
276
277 enum
278 {
279     enbnxninitcombruleDETECT,
280     enbnxninitcombruleGEOM,
281     enbnxninitcombruleLB,
282     enbnxninitcombruleNONE
283 };
284
285 /* Initialize the non-bonded atom data structure.
286  * The enum for nbatXFormat is in the file defining nbnxn_atomdata_t.
287  * Copy the ntypes*ntypes*2 sized nbfp non-bonded parameter list
288  * to the atom data structure.
289  * enbnxninitcombrule sets what combination rule data gets stored in nbat.
290  */
291 void nbnxn_atomdata_init(const gmx::MDLogger& mdlog,
292                          nbnxn_atomdata_t*    nbat,
293                          Nbnxm::KernelType    kernelType,
294                          int                  enbnxninitcombrule,
295                          int                  ntype,
296                          const real*          nbfp,
297                          int                  n_energygroups,
298                          int                  nout);
299
300 void nbnxn_atomdata_set(nbnxn_atomdata_t*     nbat,
301                         const Nbnxm::GridSet& gridSet,
302                         const t_mdatoms*      mdatoms,
303                         const int*            atinfo);
304
305 /* Copy the shift vectors to nbat */
306 void nbnxn_atomdata_copy_shiftvec(gmx_bool dynamic_box, rvec* shift_vec, nbnxn_atomdata_t* nbat);
307
308 /*! \brief Transform coordinates to xbat layout
309  *
310  * Creates a copy of the coordinates buffer using short-range ordering.
311  *
312  * \param[in] gridSet      The grids data.
313  * \param[in] locality     If the transformation should be applied to local or non local coordinates.
314  * \param[in] fillLocal    Tells if the local filler particle coordinates should be zeroed.
315  * \param[in] coordinates  Coordinates in plain rvec format.
316  * \param[in,out] nbat     Data in NBNXM format, used for mapping formats and to locate the output buffer.
317  */
318 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet& gridSet,
319                                      gmx::AtomLocality     locality,
320                                      bool                  fillLocal,
321                                      const rvec*           coordinates,
322                                      nbnxn_atomdata_t*     nbat);
323
324 /*! \brief Transform coordinates to xbat layout on GPU
325  *
326  * Creates a GPU copy of the coordinates buffer using short-range ordering.
327  * As input, uses coordinates in plain rvec format in GPU memory.
328  *
329  * \param[in]     gridSet    The grids data.
330  * \param[in]     locality   If the transformation should be applied to local or non local coordinates.
331  * \param[in]     fillLocal  Tells if the local filler particle coordinates should be zeroed.
332  * \param[in,out] gpu_nbv    The NBNXM GPU data structure.
333  * \param[in]     d_x        Coordinates to be copied (in plain rvec format).
334  * \param[in]     xReadyOnDevice   Event synchronizer indicating that the coordinates are ready in the device memory.
335  */
336 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet& gridSet,
337                                     gmx::AtomLocality     locality,
338                                     bool                  fillLocal,
339                                     gmx_nbnxn_gpu_t*      gpu_nbv,
340                                     DeviceBuffer<float>   d_x,
341                                     GpuEventSynchronizer* xReadyOnDevice);
342
343 /*! \brief Add the computed forces to \p f, an internal reduction might be performed as well
344  *
345  * \param[in]  nbat        Atom data in NBNXM format.
346  * \param[in]  locality    If the reduction should be performed on local or non-local atoms.
347  * \param[in]  gridSet     The grids data.
348  * \param[out] totalForce  Buffer to accumulate resulting force
349  */
350 void reduceForces(nbnxn_atomdata_t* nbat, gmx::AtomLocality locality, const Nbnxm::GridSet& gridSet, rvec* totalForce);
351
352 /*! \brief Reduce forces on the GPU
353  *
354  * \param[in]  locality             If the reduction should be performed on local or non-local atoms.
355  * \param[out] totalForcesDevice    Device buffer to accumulate resulting force.
356  * \param[in]  gridSet              The grids data.
357  * \param[in]  pmeForcesDevice      Device buffer with PME forces.
358  * \param[in]  dependencyList       List of synchronizers that represent the dependencies the reduction task needs to sync on.
359  * \param[in]  gpu_nbv              The NBNXM GPU data structure.
360  * \param[in]  useGpuFPmeReduction  Whether PME forces should be added.
361  * \param[in]  accumulateForce      Whether there are usefull data already in the total force buffer.
362  */
363 void reduceForcesGpu(gmx::AtomLocality                          locality,
364                      DeviceBuffer<float>                        totalForcesDevice,
365                      const Nbnxm::GridSet&                      gridSet,
366                      void*                                      pmeForcesDevice,
367                      gmx::ArrayRef<GpuEventSynchronizer* const> dependencyList,
368                      gmx_nbnxn_gpu_t*                           gpu_nbv,
369                      bool                                       useGpuFPmeReduction,
370                      bool                                       accumulateForce);
371
372 /* Add the fshift force stored in nbat to fshift */
373 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t& nbat, gmx::ArrayRef<gmx::RVec> fshift);
374
375 /* Get the atom start index and number of atoms for a given locality */
376 void nbnxn_get_atom_range(gmx::AtomLocality     atomLocality,
377                           const Nbnxm::GridSet& gridSet,
378                           int*                  atomStart,
379                           int*                  nAtoms);
380 #endif