d59e488e4428ee9028f8ff9fe08a0c19e1742493
[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,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_NBNXM_PAIRLIST_H
37 #define GMX_NBNXM_PAIRLIST_H
38
39 #include <cstddef>
40
41 #include "gromacs/gpu_utils/hostallocator.h"
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/mdtypes/locality.h"
44 #include "gromacs/mdtypes/nblist.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 // This file with constants is separate from this file to be able
51 // to include it during OpenCL jitting without including config.h
52 #include "constants.h"
53 #include "pairlistparams.h"
54
55 struct NbnxnPairlistCpuWork;
56 struct NbnxnPairlistGpuWork;
57
58
59 /* Convenience type for vector with aligned memory */
60 template<typename T>
61 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
62
63 /* Convenience type for vector that avoids initialization at resize() */
64 template<typename T>
65 using FastVector = std::vector<T, gmx::DefaultInitializationAllocator<T>>;
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     int dummy[16];
75 } gmx_cache_protect_t;
76
77 /* This is the actual cluster-pair list j-entry.
78  * cj is the j-cluster.
79  * The interaction bits in excl are indexed i-major, j-minor.
80  * The cj entries are sorted such that ones with exclusions come first.
81  * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
82  * is found, all subsequent j-entries in the i-entry also have full masks.
83  */
84 struct nbnxn_cj_t
85 {
86     int          cj;   /* The j-cluster                    */
87     unsigned int excl; /* The exclusion (interaction) bits */
88 };
89
90 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
91  * The upper bits contain information for non-bonded kernel optimization.
92  * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
93  * But three flags can be used to skip interactions, currently only for subc=0
94  * !(shift & NBNXN_CI_DO_LJ(subc))   => we can skip LJ for all pairs
95  * shift & NBNXN_CI_HALF_LJ(subc)    => we can skip LJ for the second half of i
96  * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
97  */
98 #define NBNXN_CI_SHIFT 127
99 #define NBNXN_CI_DO_LJ(subc) (1 << (7 + 3 * (subc)))
100 #define NBNXN_CI_HALF_LJ(subc) (1 << (8 + 3 * (subc)))
101 #define NBNXN_CI_DO_COUL(subc) (1 << (9 + 3 * (subc)))
102
103 /* Cluster-pair Interaction masks
104  * Bit i*j-cluster-size + j tells if atom i and j interact.
105  */
106 // TODO: Rename according to convention when moving into Nbnxn namespace
107 /* All interaction mask is the same for all kernels */
108 constexpr unsigned int NBNXN_INTERACTION_MASK_ALL = 0xffffffffU;
109 /* 4x4 kernel diagonal mask */
110 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG = 0x08ceU;
111 /* 4x2 kernel diagonal masks */
112 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_0 = 0x0002U;
113 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_1 = 0x002fU;
114 /* 4x8 kernel diagonal masks */
115 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_0 = 0xf0f8fcfeU;
116 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_1 = 0x0080c0e0U;
117
118 /* Simple pair-list i-unit */
119 struct nbnxn_ci_t
120 {
121     int ci;           /* i-cluster             */
122     int shift;        /* Shift vector index plus possible flags, see above */
123     int cj_ind_start; /* Start index into cj   */
124     int cj_ind_end;   /* End index into cj     */
125 };
126
127 /* Grouped pair-list i-unit */
128 typedef struct
129 {
130     /* Returns the number of j-cluster groups in this entry */
131     int numJClusterGroups() const { return cj4_ind_end - cj4_ind_start; }
132
133     int sci;           /* i-super-cluster       */
134     int shift;         /* Shift vector index plus possible flags */
135     int cj4_ind_start; /* Start index into cj4  */
136     int cj4_ind_end;   /* End index into cj4    */
137 } nbnxn_sci_t;
138
139 /* Interaction data for a j-group for one warp */
140 struct nbnxn_im_ei_t
141 {
142     // The i-cluster interactions mask for 1 warp
143     unsigned int imask = 0U;
144     // Index into the exclusion array for 1 warp, default index 0 which means no exclusions
145     int excl_ind = 0;
146 };
147
148 typedef struct
149 {
150     int           cj[c_nbnxnGpuJgroupSize];         /* The 4 j-clusters */
151     nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data       for 2 warps */
152 } nbnxn_cj4_t;
153
154 /* Struct for storing the atom-pair interaction bits for a cluster pair in a GPU pairlist */
155 struct nbnxn_excl_t
156 {
157     /* Constructor, sets no exclusions, so all atom pairs interacting */
158     nbnxn_excl_t()
159     {
160         for (unsigned int& pairEntry : pair)
161         {
162             pairEntry = NBNXN_INTERACTION_MASK_ALL;
163         }
164     }
165
166     /* Topology exclusion interaction bits per warp */
167     unsigned int pair[c_nbnxnGpuExclSize];
168 };
169
170 /* Cluster pairlist type for use on CPUs */
171 struct NbnxnPairlistCpu
172 {
173     NbnxnPairlistCpu();
174
175     gmx_cache_protect_t cp0;
176
177     int                    na_ci;   /* The number of atoms per i-cluster        */
178     int                    na_cj;   /* The number of atoms per j-cluster        */
179     real                   rlist;   /* The radius for constructing the list     */
180     FastVector<nbnxn_ci_t> ci;      /* The i-cluster list                       */
181     FastVector<nbnxn_ci_t> ciOuter; /* The outer, unpruned i-cluster list       */
182
183     FastVector<nbnxn_cj_t> cj;      /* The j-cluster list, size ncj             */
184     FastVector<nbnxn_cj_t> cjOuter; /* The outer, unpruned j-cluster list       */
185     int                    ncjInUse; /* The number of j-clusters that are used by ci entries in this list, will be <= cj.size() */
186
187     int nci_tot; /* The total number of i clusters           */
188
189     /* Working data storage for list construction */
190     std::unique_ptr<NbnxnPairlistCpuWork> work;
191
192     gmx_cache_protect_t cp1;
193 };
194
195 /* Cluster pairlist type, with extra hierarchies, for on the GPU
196  *
197  * NOTE: for better performance when combining lists over threads,
198  *       all vectors should use default initialization. But when
199  *       changing this, excl should be intialized when adding entries.
200  */
201 struct NbnxnPairlistGpu
202 {
203     /* Constructor
204      *
205      * \param[in] pinningPolicy  Sets the pinning policy for all buffers used on the GPU
206      */
207     NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy);
208
209     gmx_cache_protect_t cp0;
210
211     int  na_ci; /* The number of atoms per i-cluster        */
212     int  na_cj; /* The number of atoms per j-cluster        */
213     int  na_sc; /* The number of atoms per super cluster    */
214     real rlist; /* The radius for constructing the list     */
215     // The i-super-cluster list, indexes into cj4;
216     gmx::HostVector<nbnxn_sci_t> sci;
217     // The list of 4*j-cluster groups
218     gmx::HostVector<nbnxn_cj4_t> cj4;
219     // Atom interaction bits (non-exclusions)
220     gmx::HostVector<nbnxn_excl_t> excl;
221     // The total number of i-clusters
222     int nci_tot;
223
224     /* Working data storage for list construction */
225     std::unique_ptr<NbnxnPairlistGpuWork> work;
226
227     gmx_cache_protect_t cp1;
228 };
229
230 //! Initializes a free-energy pair-list
231 void nbnxn_init_pairlist_fep(t_nblist* nl);
232
233 #endif