SYCL: Use acc.bind(cgh) instead of cgh.require(acc)
[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 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
37 #ifndef GMX_NBNXM_PAIRLIST_H
38 #define GMX_NBNXM_PAIRLIST_H
39
40 #include <cstddef>
41
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/defaultinitializationallocator.h"
47 #include "gromacs/utility/enumerationhelpers.h"
48 #include "gromacs/utility/real.h"
49
50 #include "pairlistparams.h"
51
52 struct NbnxnPairlistCpuWork;
53 struct NbnxnPairlistGpuWork;
54 struct t_nblist;
55
56
57 //! Convenience type for vector with aligned memory
58 template<typename T>
59 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
60
61 //! Convenience type for vector that avoids initialization at resize()
62 template<typename T>
63 using FastVector = std::vector<T, gmx::DefaultInitializationAllocator<T>>;
64
65 /*! \brief Cache-line protection buffer
66  *
67  * A buffer data structure of 64 bytes
68  * to be placed at the beginning and end of structs
69  * to avoid cache invalidation of the real contents
70  * of the struct by writes to neighboring memory.
71  */
72 typedef struct
73 {
74     //! Unused field used to create space to protect cache lines that are in use
75     int dummy[16];
76 } gmx_cache_protect_t;
77
78 /*! \brief This is the actual cluster-pair list j-entry.
79  *
80  * cj is the j-cluster.
81  * The interaction bits in excl are indexed i-major, j-minor.
82  * The cj entries are sorted such that ones with exclusions come first.
83  * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
84  * is found, all subsequent j-entries in the i-entry also have full masks.
85  */
86 struct nbnxn_cj_t
87 {
88     //! The j-cluster
89     int cj;
90     //! The exclusion (interaction) bits
91     unsigned int excl;
92 };
93
94 /*! \brief Constants for interpreting interaction flags
95  *
96  * In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
97  * The upper bits contain information for non-bonded kernel optimization.
98  * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
99  * But three flags can be used to skip interactions, currently only for subc=0
100  * !(shift & NBNXN_CI_DO_LJ(subc))   => we can skip LJ for all pairs
101  * shift & NBNXN_CI_HALF_LJ(subc)    => we can skip LJ for the second half of i
102  * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
103  */
104 //! \{
105 #define NBNXN_CI_SHIFT 127
106 #define NBNXN_CI_DO_LJ(subc) (1 << (7 + 3 * (subc)))
107 #define NBNXN_CI_HALF_LJ(subc) (1 << (8 + 3 * (subc)))
108 #define NBNXN_CI_DO_COUL(subc) (1 << (9 + 3 * (subc)))
109 //! \}
110
111 /*! \brief Cluster-pair Interaction masks
112  *
113  * Bit i*j-cluster-size + j tells if atom i and j interact.
114  */
115 //! \{
116 // TODO: Rename according to convention when moving into Nbnxn namespace
117 //! All interaction mask is the same for all kernels
118 constexpr unsigned int NBNXN_INTERACTION_MASK_ALL = 0xffffffffU;
119 //! 4x4 kernel diagonal mask
120 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG = 0x08ceU;
121 //! 4x2 kernel diagonal masks
122 //! \{
123 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_0 = 0x0002U;
124 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_1 = 0x002fU;
125 //! \}
126 //! 4x8 kernel diagonal masks
127 //! \{
128 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_0 = 0xf0f8fcfeU;
129 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_1 = 0x0080c0e0U;
130 //! \}
131 //! \}
132
133 /*! \brief Lower limit for square interaction distances in nonbonded kernels.
134  *
135  * For smaller values we will overflow when calculating r^-1 or r^-12, but
136  * to keep it simple we always apply the limit from the tougher r^-12 condition.
137  */
138 #if GMX_DOUBLE
139 // Some double precision SIMD architectures use single precision in the first
140 // step, so although the double precision criterion would allow smaller rsq,
141 // we need to stay in single precision with some margin for the N-R iterations.
142 constexpr double c_nbnxnMinDistanceSquared = 1.0e-36;
143 #else
144 // The worst intermediate value we might evaluate is r^-12, which
145 // means we should ensure r^2 stays above pow(GMX_FLOAT_MAX,-1.0/6.0)*1.01 (some margin)
146 constexpr float c_nbnxnMinDistanceSquared = 3.82e-07F; // r > 6.2e-4
147 #endif
148
149
150 //! The number of clusters in a super-cluster, used for GPU
151 constexpr int c_nbnxnGpuNumClusterPerSupercluster = 8;
152
153 /*! \brief With GPU kernels we group cluster pairs in 4 to optimize memory usage
154  * of integers containing 32 bits.
155  */
156 constexpr int c_nbnxnGpuJgroupSize = (32 / c_nbnxnGpuNumClusterPerSupercluster);
157
158 /*! \internal
159  * \brief Simple pair-list i-unit
160  */
161 struct nbnxn_ci_t
162 {
163     //! i-cluster
164     int ci;
165     //! Shift vector index plus possible flags, see above
166     int shift;
167     //! Start index into cj
168     int cj_ind_start;
169     //! End index into cj
170     int cj_ind_end;
171 };
172
173 //! Grouped pair-list i-unit
174 typedef struct nbnxn_sci
175 {
176     //! Returns the number of j-cluster groups in this entry
177     int numJClusterGroups() const { return cj4_ind_end - cj4_ind_start; }
178
179     //! i-super-cluster
180     int sci;
181     //! Shift vector index plus possible flags
182     int shift;
183     //! Start index into cj4
184     int cj4_ind_start;
185     //! End index into cj4
186     int cj4_ind_end;
187 } nbnxn_sci_t;
188
189 //! Interaction data for a j-group for one warp
190 struct nbnxn_im_ei_t
191 {
192     //! The i-cluster interactions mask for 1 warp
193     unsigned int imask = 0U;
194     //! Index into the exclusion array for 1 warp, default index 0 which means no exclusions
195     int excl_ind = 0;
196 };
197
198 //! Four-way j-cluster lists
199 typedef struct
200 {
201     //! The 4 j-clusters
202     int cj[c_nbnxnGpuJgroupSize];
203     //! The i-cluster mask data for 2 warps
204     nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit];
205 } nbnxn_cj4_t;
206
207 //! Struct for storing the atom-pair interaction bits for a cluster pair in a GPU pairlist
208 struct nbnxn_excl_t
209 {
210     //! Constructor, sets no exclusions, so all atom pairs interacting
211     MSVC_DIAGNOSTIC_IGNORE(26495) // pair is not being initialized!
212     nbnxn_excl_t()
213     {
214         for (unsigned int& pairEntry : pair)
215         {
216             pairEntry = NBNXN_INTERACTION_MASK_ALL;
217         }
218     }
219     MSVC_DIAGNOSTIC_RESET
220
221     //! Topology exclusion interaction bits per warp
222     unsigned int pair[c_nbnxnGpuExclSize];
223 };
224
225 //! Cluster pairlist type for use on CPUs
226 struct NbnxnPairlistCpu
227 {
228     NbnxnPairlistCpu();
229
230     //! Cache protection
231     gmx_cache_protect_t cp0;
232
233     //! The number of atoms per i-cluster
234     int na_ci;
235     //! The number of atoms per j-cluster
236     int na_cj;
237     //! The radius for constructing the list
238     real rlist;
239     //! The i-cluster list
240     FastVector<nbnxn_ci_t> ci;
241     //! The outer, unpruned i-cluster list
242     FastVector<nbnxn_ci_t> ciOuter;
243
244     //! The j-cluster list, size ncj
245     FastVector<nbnxn_cj_t> cj;
246     //! The outer, unpruned j-cluster list
247     FastVector<nbnxn_cj_t> cjOuter;
248     //! The number of j-clusters that are used by ci entries in this list, will be <= cj.size()
249     int ncjInUse;
250
251     //! The total number of i clusters
252     int nci_tot;
253
254     //! Working data storage for list construction
255     std::unique_ptr<NbnxnPairlistCpuWork> work;
256
257     //! Cache protection
258     gmx_cache_protect_t cp1;
259 };
260
261 /* Cluster pairlist type, with extra hierarchies, for on the GPU
262  *
263  * NOTE: for better performance when combining lists over threads,
264  *       all vectors should use default initialization. But when
265  *       changing this, excl should be intialized when adding entries.
266  */
267 struct NbnxnPairlistGpu
268 {
269     /*! \brief Constructor
270      *
271      * \param[in] pinningPolicy  Sets the pinning policy for all buffers used on the GPU
272      */
273     NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy);
274
275     //! Cache protection
276     gmx_cache_protect_t cp0;
277
278     //! The number of atoms per i-cluster
279     int na_ci;
280     //! The number of atoms per j-cluster
281     int na_cj;
282     //! The number of atoms per super cluster
283     int na_sc;
284     //! The radius for constructing the list
285     real rlist;
286     // The i-super-cluster list, indexes into cj4;
287     gmx::HostVector<nbnxn_sci_t> sci;
288     // The list of 4*j-cluster groups
289     gmx::HostVector<nbnxn_cj4_t> cj4;
290     // Atom interaction bits (non-exclusions)
291     gmx::HostVector<nbnxn_excl_t> excl;
292     // The total number of i-clusters
293     int nci_tot;
294
295     //! Working data storage for list construction
296     std::unique_ptr<NbnxnPairlistGpuWork> work;
297
298     //! Cache protection
299     gmx_cache_protect_t cp1;
300 };
301
302 #endif