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