Split lines with many copyright years
[gromacs.git] / src / gromacs / nbnxm / pairlist.h
blob777a42a5a0c4c55727202a8276a63b39e65f86a8
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020, 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.
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.
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.
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.
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.
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.
37 #ifndef GMX_NBNXM_PAIRLIST_H
38 #define GMX_NBNXM_PAIRLIST_H
40 #include <cstddef>
42 #include "gromacs/gpu_utils/hostallocator.h"
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/mdtypes/locality.h"
45 #include "gromacs/mdtypes/nblist.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/defaultinitializationallocator.h"
48 #include "gromacs/utility/enumerationhelpers.h"
49 #include "gromacs/utility/real.h"
51 // This file with constants is separate from this file to be able
52 // to include it during OpenCL jitting without including config.h
53 #include "constants.h"
54 #include "pairlistparams.h"
56 struct NbnxnPairlistCpuWork;
57 struct NbnxnPairlistGpuWork;
60 /* Convenience type for vector with aligned memory */
61 template<typename T>
62 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
64 /* Convenience type for vector that avoids initialization at resize() */
65 template<typename T>
66 using FastVector = std::vector<T, gmx::DefaultInitializationAllocator<T>>;
68 /* A buffer data structure of 64 bytes
69 * to be placed at the beginning and end of structs
70 * to avoid cache invalidation of the real contents
71 * of the struct by writes to neighboring memory.
73 typedef struct
75 int dummy[16];
76 } gmx_cache_protect_t;
78 /* This is the actual cluster-pair list j-entry.
79 * cj is the j-cluster.
80 * The interaction bits in excl are indexed i-major, j-minor.
81 * The cj entries are sorted such that ones with exclusions come first.
82 * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
83 * is found, all subsequent j-entries in the i-entry also have full masks.
85 struct nbnxn_cj_t
87 int cj; /* The j-cluster */
88 unsigned int excl; /* The exclusion (interaction) bits */
91 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
92 * The upper bits contain information for non-bonded kernel optimization.
93 * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
94 * But three flags can be used to skip interactions, currently only for subc=0
95 * !(shift & NBNXN_CI_DO_LJ(subc)) => we can skip LJ for all pairs
96 * shift & NBNXN_CI_HALF_LJ(subc) => we can skip LJ for the second half of i
97 * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
99 #define NBNXN_CI_SHIFT 127
100 #define NBNXN_CI_DO_LJ(subc) (1 << (7 + 3 * (subc)))
101 #define NBNXN_CI_HALF_LJ(subc) (1 << (8 + 3 * (subc)))
102 #define NBNXN_CI_DO_COUL(subc) (1 << (9 + 3 * (subc)))
104 /* Cluster-pair Interaction masks
105 * Bit i*j-cluster-size + j tells if atom i and j interact.
107 // TODO: Rename according to convention when moving into Nbnxn namespace
108 /* All interaction mask is the same for all kernels */
109 constexpr unsigned int NBNXN_INTERACTION_MASK_ALL = 0xffffffffU;
110 /* 4x4 kernel diagonal mask */
111 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG = 0x08ceU;
112 /* 4x2 kernel diagonal masks */
113 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_0 = 0x0002U;
114 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_1 = 0x002fU;
115 /* 4x8 kernel diagonal masks */
116 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_0 = 0xf0f8fcfeU;
117 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_1 = 0x0080c0e0U;
119 /* Simple pair-list i-unit */
120 struct nbnxn_ci_t
122 int ci; /* i-cluster */
123 int shift; /* Shift vector index plus possible flags, see above */
124 int cj_ind_start; /* Start index into cj */
125 int cj_ind_end; /* End index into cj */
128 /* Grouped pair-list i-unit */
129 typedef struct
131 /* Returns the number of j-cluster groups in this entry */
132 int numJClusterGroups() const { return cj4_ind_end - cj4_ind_start; }
134 int sci; /* i-super-cluster */
135 int shift; /* Shift vector index plus possible flags */
136 int cj4_ind_start; /* Start index into cj4 */
137 int cj4_ind_end; /* End index into cj4 */
138 } nbnxn_sci_t;
140 /* Interaction data for a j-group for one warp */
141 struct nbnxn_im_ei_t
143 // The i-cluster interactions mask for 1 warp
144 unsigned int imask = 0U;
145 // Index into the exclusion array for 1 warp, default index 0 which means no exclusions
146 int excl_ind = 0;
149 typedef struct
151 int cj[c_nbnxnGpuJgroupSize]; /* The 4 j-clusters */
152 nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data for 2 warps */
153 } nbnxn_cj4_t;
155 /* Struct for storing the atom-pair interaction bits for a cluster pair in a GPU pairlist */
156 struct nbnxn_excl_t
158 /* Constructor, sets no exclusions, so all atom pairs interacting */
159 nbnxn_excl_t()
161 for (unsigned int& pairEntry : pair)
163 pairEntry = NBNXN_INTERACTION_MASK_ALL;
167 /* Topology exclusion interaction bits per warp */
168 unsigned int pair[c_nbnxnGpuExclSize];
171 /* Cluster pairlist type for use on CPUs */
172 struct NbnxnPairlistCpu
174 NbnxnPairlistCpu();
176 gmx_cache_protect_t cp0;
178 int na_ci; /* The number of atoms per i-cluster */
179 int na_cj; /* The number of atoms per j-cluster */
180 real rlist; /* The radius for constructing the list */
181 FastVector<nbnxn_ci_t> ci; /* The i-cluster list */
182 FastVector<nbnxn_ci_t> ciOuter; /* The outer, unpruned i-cluster list */
184 FastVector<nbnxn_cj_t> cj; /* The j-cluster list, size ncj */
185 FastVector<nbnxn_cj_t> cjOuter; /* The outer, unpruned j-cluster list */
186 int ncjInUse; /* The number of j-clusters that are used by ci entries in this list, will be <= cj.size() */
188 int nci_tot; /* The total number of i clusters */
190 /* Working data storage for list construction */
191 std::unique_ptr<NbnxnPairlistCpuWork> work;
193 gmx_cache_protect_t cp1;
196 /* Cluster pairlist type, with extra hierarchies, for on the GPU
198 * NOTE: for better performance when combining lists over threads,
199 * all vectors should use default initialization. But when
200 * changing this, excl should be intialized when adding entries.
202 struct NbnxnPairlistGpu
204 /* Constructor
206 * \param[in] pinningPolicy Sets the pinning policy for all buffers used on the GPU
208 NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy);
210 gmx_cache_protect_t cp0;
212 int na_ci; /* The number of atoms per i-cluster */
213 int na_cj; /* The number of atoms per j-cluster */
214 int na_sc; /* The number of atoms per super cluster */
215 real rlist; /* The radius for constructing the list */
216 // The i-super-cluster list, indexes into cj4;
217 gmx::HostVector<nbnxn_sci_t> sci;
218 // The list of 4*j-cluster groups
219 gmx::HostVector<nbnxn_cj4_t> cj4;
220 // Atom interaction bits (non-exclusions)
221 gmx::HostVector<nbnxn_excl_t> excl;
222 // The total number of i-clusters
223 int nci_tot;
225 /* Working data storage for list construction */
226 std::unique_ptr<NbnxnPairlistGpuWork> work;
228 gmx_cache_protect_t cp1;
231 //! Initializes a free-energy pair-list
232 void nbnxn_init_pairlist_fep(t_nblist* nl);
234 #endif