Remove tree-reduce algorithm
[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 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \libinternal \file
37  *  \brief
38  *  Functionality for per-atom data in the nbnxm module
39  *
40  *  \author Berk Hess <hess@kth.se>
41  *  \ingroup module_nbnxm
42  *  \inlibraryapi
43  */
44
45
46 #ifndef GMX_NBNXN_ATOMDATA_H
47 #define GMX_NBNXN_ATOMDATA_H
48
49 #include <cstdio>
50
51 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
52 #include "gromacs/gpu_utils/hostallocator.h"
53 #include "gromacs/math/vectypes.h"
54 #include "gromacs/mdtypes/locality.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/bitmask.h"
57 #include "gromacs/utility/real.h"
58
59 namespace gmx
60 {
61 class MDLogger;
62 }
63
64 struct NbnxmGpu;
65 struct nbnxn_atomdata_t;
66 struct nonbonded_verlet_t;
67
68 class GpuEventSynchronizer;
69
70 namespace Nbnxm
71 {
72 class GridSet;
73 enum class KernelType;
74 } // namespace Nbnxm
75
76 //! Convenience type for vector with aligned memory
77 template<typename T>
78 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
79
80 enum
81 {
82     nbatXYZ,
83     nbatXYZQ,
84     nbatX4,
85     nbatX8
86 };
87
88 //! Stride for coordinate/force arrays with xyz coordinate storage
89 static constexpr int STRIDE_XYZ = 3;
90 //! Stride for coordinate/force arrays with xyzq coordinate storage
91 static constexpr int STRIDE_XYZQ = 4;
92 //! Size of packs of x, y or z with SIMD 4-grouped packed coordinates/forces
93 static constexpr int c_packX4 = 4;
94 //! Size of packs of x, y or z with SIMD 8-grouped packed coordinates/forces
95 static constexpr int c_packX8 = 8;
96 //! Stridefor a pack of 4 coordinates/forces
97 static constexpr int STRIDE_P4 = DIM * c_packX4;
98 //! Stridefor a pack of 8 coordinates/forces
99 static constexpr int STRIDE_P8 = DIM * c_packX8;
100
101 //! Returns the index in a coordinate array corresponding to atom a
102 template<int packSize>
103 static inline int atom_to_x_index(int a)
104 {
105     return DIM * (a & ~(packSize - 1)) + (a & (packSize - 1));
106 }
107
108 /*! \internal
109  * \brief Struct that holds force and energy output buffers */
110 struct nbnxn_atomdata_output_t
111 {
112     /*! \brief Constructor
113      *
114      * \param[in] kernelType              Type of non-bonded kernel
115      * \param[in] numEnergyGroups         The number of energy groups
116      * \param[in] simdEnergyBufferStride  Stride for entries in the energy buffers for SIMD kernels
117      * \param[in] pinningPolicy           Sets the pinning policy for all buffers used on the GPU
118      */
119     nbnxn_atomdata_output_t(Nbnxm::KernelType  kernelType,
120                             int                numEnergyGroups,
121                             int                simdEnergyBufferStride,
122                             gmx::PinningPolicy pinningPolicy);
123
124     //! f, size natoms*fstride
125     gmx::HostVector<real> f;
126     //! Shift force array, size SHIFTS*DIM
127     gmx::HostVector<real> fshift;
128     //! Temporary Van der Waals group energy storage
129     gmx::HostVector<real> Vvdw;
130     //! Temporary Coulomb group energy storage
131     gmx::HostVector<real> Vc;
132     //! Temporary SIMD Van der Waals group energy storage
133     AlignedVector<real> VSvdw;
134     //! Temporary SIMD Coulomb group energy storage
135     AlignedVector<real> VSc;
136 };
137
138 /*! \brief Block size in atoms for the non-bonded thread force-buffer reduction.
139  *
140  * Should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
141  * Should be small to reduce the reduction and zeroing cost,
142  * but too small will result in overhead.
143  * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
144  */
145 #if GMX_DOUBLE
146 #    define NBNXN_BUFFERFLAG_SIZE 8
147 #else
148 #    define NBNXN_BUFFERFLAG_SIZE 16
149 #endif
150
151 /*! \brief We store the reduction flags as gmx_bitmask_t.
152  * This limits the number of flags to BITMASK_SIZE.
153  */
154 #define NBNXN_BUFFERFLAG_MAX_THREADS (BITMASK_SIZE)
155
156
157 //! LJ combination rules
158 enum class LJCombinationRule : int
159 {
160     //! Geometric
161     Geometric,
162     //! Lorentz-Berthelot
163     LorentzBerthelot,
164     //! No rule
165     None,
166     //! Size of the enum
167     Count
168 };
169
170 //! String corresponding to LJ combination rule
171 const char* enumValueToString(LJCombinationRule enumValue);
172
173 /*! \internal
174  * \brief Struct that stores atom related data for the nbnxn module
175  *
176  * Note: performance would improve slightly when all std::vector containers
177  *       in this struct would not initialize during resize().
178  */
179 struct nbnxn_atomdata_t
180 { //NOLINT(clang-analyzer-optin.performance.Padding)
181     /*! \internal
182      * \brief The actual atom data parameter values */
183     struct Params
184     {
185         /*! \brief Constructor
186          *
187          * \param[in] pinningPolicy  Sets the pinning policy for all data that might be transfered to a GPU
188          */
189         Params(gmx::PinningPolicy pinningPolicy);
190
191         //! The number of different atom types
192         int numTypes;
193         //! Lennard-Jone 6*C6 and 12*C12 parameters, size numTypes*2*2
194         gmx::HostVector<real> nbfp;
195         //! Combination rule, see enum defined above
196         LJCombinationRule ljCombinationRule;
197         //! LJ parameters per atom type, size numTypes*2
198         gmx::HostVector<real> nbfp_comb;
199         //! As nbfp, but with a stride for the present SIMD architecture
200         AlignedVector<real> nbfp_aligned;
201         //! Atom types per atom
202         gmx::HostVector<int> type;
203         //! LJ parameters per atom for fast SIMD loading
204         gmx::HostVector<real> lj_comb;
205         //! Charges per atom, not set with format nbatXYZQ
206         gmx::HostVector<real> q;
207         //! The number of energy groups
208         int nenergrp;
209         //! 2log(nenergrp)
210         int neg_2log;
211         //! The energy groups, one int entry per cluster, only set when needed
212         gmx::HostVector<int> energrp;
213     };
214
215     /*! \internal
216      * \brief Diagonal and topology exclusion helper data for all SIMD kernels. */
217     struct SimdMasks
218     {
219         SimdMasks();
220
221         //! Helper data for setting up diagonal exclusion masks in the SIMD 4xN kernels
222         AlignedVector<real> diagonal_4xn_j_minus_i;
223         //! Helper data for setting up diaginal exclusion masks in the SIMD 2xNN kernels
224         AlignedVector<real> diagonal_2xnn_j_minus_i;
225         //! Filters for topology exclusion masks for the SIMD kernels
226         AlignedVector<uint32_t> exclusion_filter;
227         //! Filters for topology exclusion masks for double SIMD kernels without SIMD int32 logical support
228         AlignedVector<uint64_t> exclusion_filter64;
229         //! Array of masks needed for exclusions
230         AlignedVector<real> interaction_array;
231     };
232
233     /*! \brief Constructor
234      *
235      * \param[in] pinningPolicy  Sets the pinning policy for all data that might be transfered to a GPU
236      */
237     nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy);
238
239     //! Returns a const reference to the parameters
240     const Params& params() const { return params_; }
241
242     //! Returns a non-const reference to the parameters
243     Params& paramsDeprecated() { return params_; }
244
245     //! Returns the current total number of atoms stored
246     int numAtoms() const { return numAtoms_; }
247
248     //! Return the coordinate buffer, and q with xFormat==nbatXYZQ
249     gmx::ArrayRef<const real> x() const { return x_; }
250
251     //! Return the coordinate buffer, and q with xFormat==nbatXYZQ
252     gmx::ArrayRef<real> x() { return x_; }
253
254     //! Resizes the coordinate buffer and sets the number of atoms
255     void resizeCoordinateBuffer(int numAtoms);
256
257     //! Resizes the force buffers for the current number of atoms
258     void resizeForceBuffers();
259
260 private:
261     //! The LJ and charge parameters
262     Params params_;
263     //! The total number of atoms currently stored
264     int numAtoms_;
265
266 public:
267     //! Number of local atoms
268     int natoms_local;
269     //! The format of x (and q), enum
270     int XFormat;
271     //! The format of f, enum
272     int FFormat;
273     //! Do we need to update shift_vec every step?
274     gmx_bool bDynamicBox;
275     //! Shift vectors, copied from t_forcerec
276     gmx::HostVector<gmx::RVec> shift_vec;
277     //! stride for a coordinate in x (usually 3 or 4)
278     int xstride;
279     //! stride for a coordinate in f (usually 3 or 4)
280     int fstride;
281
282 private:
283     //! x and possibly q, size natoms*xstride
284     gmx::HostVector<real> x_;
285
286 public:
287     //! Masks for handling exclusions in the SIMD kernels
288     const SimdMasks simdMasks;
289
290     //! Output data structures, 1 per thread
291     std::vector<nbnxn_atomdata_output_t> out;
292
293     //! Reduction related data
294     //! \{
295     //! Use the flags or operate on all atoms
296     gmx_bool bUseBufferFlags;
297     //! Flags for buffer zeroing+reduc.
298     std::vector<gmx_bitmask_t> buffer_flags;
299     //! \}
300 };
301
302 /*! \brief Copy na rvec elements from x to xnb using nbatFormat, start dest a0,
303  * and fills up to na_round with coordinates that are far away.
304  */
305 void copy_rvec_to_nbat_real(const int* a, int na, int na_round, const rvec* x, int nbatFormat, real* xnb, int a0);
306
307 //! Describes the combination rule in use by this force field
308 enum
309 {
310     enbnxninitcombruleDETECT,
311     enbnxninitcombruleGEOM,
312     enbnxninitcombruleLB,
313     enbnxninitcombruleNONE
314 };
315
316 /*! \brief Initialize the non-bonded atom data structure.
317  *
318  * The enum for nbatXFormat is in the file defining nbnxn_atomdata_t.
319  * Copy the ntypes*ntypes*2 sized nbfp non-bonded parameter list
320  * to the atom data structure.
321  * enbnxninitcombrule sets what combination rule data gets stored in nbat.
322  */
323 void nbnxn_atomdata_init(const gmx::MDLogger&      mdlog,
324                          nbnxn_atomdata_t*         nbat,
325                          Nbnxm::KernelType         kernelType,
326                          int                       enbnxninitcombrule,
327                          int                       ntype,
328                          gmx::ArrayRef<const real> nbfp,
329                          int                       n_energygroups,
330                          int                       nout);
331
332 //! Sets the atomdata after pair search
333 void nbnxn_atomdata_set(nbnxn_atomdata_t*         nbat,
334                         const Nbnxm::GridSet&     gridSet,
335                         gmx::ArrayRef<const int>  atomTypes,
336                         gmx::ArrayRef<const real> atomCharges,
337                         gmx::ArrayRef<const int>  atomInfo);
338
339 //! Copy the shift vectors to nbat
340 void nbnxn_atomdata_copy_shiftvec(gmx_bool dynamic_box, rvec* shift_vec, nbnxn_atomdata_t* nbat);
341
342 /*! \brief Transform coordinates to xbat layout
343  *
344  * Creates a copy of the coordinates buffer using short-range ordering.
345  *
346  * \param[in] gridSet      The grids data.
347  * \param[in] locality     If the transformation should be applied to local or non local coordinates.
348  * \param[in] fillLocal    Tells if the local filler particle coordinates should be zeroed.
349  * \param[in] coordinates  Coordinates in plain rvec format.
350  * \param[in,out] nbat     Data in NBNXM format, used for mapping formats and to locate the output buffer.
351  */
352 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet& gridSet,
353                                      gmx::AtomLocality     locality,
354                                      bool                  fillLocal,
355                                      const rvec*           coordinates,
356                                      nbnxn_atomdata_t*     nbat);
357
358 /*! \brief Transform coordinates to xbat layout on GPU
359  *
360  * Creates a GPU copy of the coordinates buffer using short-range ordering.
361  * As input, uses coordinates in plain rvec format in GPU memory.
362  *
363  * \param[in]     gridSet    The grids data.
364  * \param[in]     locality   If the transformation should be applied to local or non local coordinates.
365  * \param[in]     fillLocal  Tells if the local filler particle coordinates should be zeroed.
366  * \param[in,out] gpu_nbv    The NBNXM GPU data structure.
367  * \param[in]     d_x        Coordinates to be copied (in plain rvec format).
368  * \param[in]     xReadyOnDevice   Event synchronizer indicating that the coordinates are ready in the device memory.
369  */
370 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet&   gridSet,
371                                     gmx::AtomLocality       locality,
372                                     bool                    fillLocal,
373                                     NbnxmGpu*               gpu_nbv,
374                                     DeviceBuffer<gmx::RVec> d_x,
375                                     GpuEventSynchronizer*   xReadyOnDevice);
376
377 /*! \brief Add the computed forces to \p f, an internal reduction might be performed as well
378  *
379  * \param[in]  nbat        Atom data in NBNXM format.
380  * \param[in]  locality    If the reduction should be performed on local or non-local atoms.
381  * \param[in]  gridSet     The grids data.
382  * \param[out] totalForce  Buffer to accumulate resulting force
383  */
384 void reduceForces(nbnxn_atomdata_t* nbat, gmx::AtomLocality locality, const Nbnxm::GridSet& gridSet, rvec* totalForce);
385
386 //! Add the fshift force stored in nbat to fshift
387 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t& nbat, gmx::ArrayRef<gmx::RVec> fshift);
388
389 //! Get the atom start index and number of atoms for a given locality
390 void nbnxn_get_atom_range(gmx::AtomLocality     atomLocality,
391                           const Nbnxm::GridSet& gridSet,
392                           int*                  atomStart,
393                           int*                  nAtoms);
394 #endif