Make eelType and evdwType scoped enums, + cleanup
[alexxy/gromacs.git] / src / gromacs / nbnxm / cuda / nbnxm_cuda.cu
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, 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 /*! \file
37  *  \brief Define CUDA implementation of nbnxn_gpu.h
38  *
39  *  \author Szilard Pall <pall.szilard@gmail.com>
40  */
41 #include "gmxpre.h"
42
43 #include "config.h"
44
45 #include <assert.h>
46 #include <stdlib.h>
47
48 #include "gromacs/nbnxm/nbnxm_gpu.h"
49
50 #if defined(_MSVC)
51 #    include <limits>
52 #endif
53
54
55 #include "nbnxm_cuda.h"
56
57 #include "gromacs/gpu_utils/gpu_utils.h"
58 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
59 #include "gromacs/gpu_utils/typecasts.cuh"
60 #include "gromacs/gpu_utils/vectype_ops.cuh"
61 #include "gromacs/hardware/device_information.h"
62 #include "gromacs/mdtypes/simulation_workload.h"
63 #include "gromacs/nbnxm/atomdata.h"
64 #include "gromacs/nbnxm/gpu_common.h"
65 #include "gromacs/nbnxm/gpu_common_utils.h"
66 #include "gromacs/nbnxm/gpu_data_mgmt.h"
67 #include "gromacs/nbnxm/grid.h"
68 #include "gromacs/nbnxm/nbnxm.h"
69 #include "gromacs/nbnxm/pairlist.h"
70 #include "gromacs/timing/gpu_timing.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/gmxassert.h"
73
74 #include "nbnxm_buffer_ops_kernels.cuh"
75 #include "nbnxm_cuda_types.h"
76
77 /***** The kernel declarations/definitions come here *****/
78
79 /* Top-level kernel declaration generation: will generate through multiple
80  * inclusion the following flavors for all kernel declarations:
81  * - force-only output;
82  * - force and energy output;
83  * - force-only with pair list pruning;
84  * - force and energy output with pair list pruning.
85  */
86 #define FUNCTION_DECLARATION_ONLY
87 /** Force only **/
88 #include "nbnxm_cuda_kernels.cuh"
89 /** Force & energy **/
90 #define CALC_ENERGIES
91 #include "nbnxm_cuda_kernels.cuh"
92 #undef CALC_ENERGIES
93
94 /*** Pair-list pruning kernels ***/
95 /** Force only **/
96 #define PRUNE_NBL
97 #include "nbnxm_cuda_kernels.cuh"
98 /** Force & energy **/
99 #define CALC_ENERGIES
100 #include "nbnxm_cuda_kernels.cuh"
101 #undef CALC_ENERGIES
102 #undef PRUNE_NBL
103
104 /* Prune-only kernels */
105 #include "nbnxm_cuda_kernel_pruneonly.cuh"
106 #undef FUNCTION_DECLARATION_ONLY
107
108 /* Now generate the function definitions if we are using a single compilation unit. */
109 #if GMX_CUDA_NB_SINGLE_COMPILATION_UNIT
110 #    include "nbnxm_cuda_kernel_F_noprune.cu"
111 #    include "nbnxm_cuda_kernel_F_prune.cu"
112 #    include "nbnxm_cuda_kernel_VF_noprune.cu"
113 #    include "nbnxm_cuda_kernel_VF_prune.cu"
114 #    include "nbnxm_cuda_kernel_pruneonly.cu"
115 #endif /* GMX_CUDA_NB_SINGLE_COMPILATION_UNIT */
116
117 namespace Nbnxm
118 {
119
120 //! Number of CUDA threads in a block
121 // TODO Optimize this through experimentation
122 constexpr static int c_bufOpsThreadsPerBlock = 128;
123
124 /*! Nonbonded kernel function pointer type */
125 typedef void (*nbnxn_cu_kfunc_ptr_t)(const cu_atomdata_t, const NBParamGpu, const gpu_plist, bool);
126
127 /*********************************/
128
129 /*! Returns the number of blocks to be used for the nonbonded GPU kernel. */
130 static inline int calc_nb_kernel_nblock(int nwork_units, const DeviceInformation* deviceInfo)
131 {
132     int max_grid_x_size;
133
134     assert(deviceInfo);
135     /* CUDA does not accept grid dimension of 0 (which can happen e.g. with an
136        empty domain) and that case should be handled before this point. */
137     assert(nwork_units > 0);
138
139     max_grid_x_size = deviceInfo->prop.maxGridSize[0];
140
141     /* do we exceed the grid x dimension limit? */
142     if (nwork_units > max_grid_x_size)
143     {
144         gmx_fatal(FARGS,
145                   "Watch out, the input system is too large to simulate!\n"
146                   "The number of nonbonded work units (=number of super-clusters) exceeds the"
147                   "maximum grid size in x dimension (%d > %d)!",
148                   nwork_units, max_grid_x_size);
149     }
150
151     return nwork_units;
152 }
153
154
155 /* Constant arrays listing all kernel function pointers and enabling selection
156    of a kernel in an elegant manner. */
157
158 /*! Pointers to the non-bonded kernels organized in 2-dim arrays by:
159  *  electrostatics and VDW type.
160  *
161  *  Note that the row- and column-order of function pointers has to match the
162  *  order of corresponding enumerated electrostatics and vdw types, resp.,
163  *  defined in nbnxn_cuda_types.h.
164  */
165
166 /*! Force-only kernel function pointers. */
167 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_noprune_ptr[c_numElecTypes][c_numVdwTypes] = {
168     { nbnxn_kernel_ElecCut_VdwLJ_F_cuda, nbnxn_kernel_ElecCut_VdwLJCombGeom_F_cuda,
169       nbnxn_kernel_ElecCut_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecCut_VdwLJFsw_F_cuda,
170       nbnxn_kernel_ElecCut_VdwLJPsw_F_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombGeom_F_cuda,
171       nbnxn_kernel_ElecCut_VdwLJEwCombLB_F_cuda },
172     { nbnxn_kernel_ElecRF_VdwLJ_F_cuda, nbnxn_kernel_ElecRF_VdwLJCombGeom_F_cuda,
173       nbnxn_kernel_ElecRF_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecRF_VdwLJFsw_F_cuda,
174       nbnxn_kernel_ElecRF_VdwLJPsw_F_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_cuda,
175       nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_cuda },
176     { nbnxn_kernel_ElecEwQSTab_VdwLJ_F_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_F_cuda,
177       nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_cuda,
178       nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_F_cuda,
179       nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_F_cuda },
180     { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_F_cuda,
181       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_cuda,
182       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_F_cuda,
183       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_F_cuda },
184     { nbnxn_kernel_ElecEw_VdwLJ_F_cuda, nbnxn_kernel_ElecEw_VdwLJCombGeom_F_cuda,
185       nbnxn_kernel_ElecEw_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecEw_VdwLJFsw_F_cuda,
186       nbnxn_kernel_ElecEw_VdwLJPsw_F_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_cuda,
187       nbnxn_kernel_ElecEw_VdwLJEwCombLB_F_cuda },
188     { nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_cuda,
189       nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_cuda,
190       nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_cuda,
191       nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_F_cuda }
192 };
193
194 /*! Force + energy kernel function pointers. */
195 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_noprune_ptr[c_numElecTypes][c_numVdwTypes] = {
196     { nbnxn_kernel_ElecCut_VdwLJ_VF_cuda, nbnxn_kernel_ElecCut_VdwLJCombGeom_VF_cuda,
197       nbnxn_kernel_ElecCut_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecCut_VdwLJFsw_VF_cuda,
198       nbnxn_kernel_ElecCut_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombGeom_VF_cuda,
199       nbnxn_kernel_ElecCut_VdwLJEwCombLB_VF_cuda },
200     { nbnxn_kernel_ElecRF_VdwLJ_VF_cuda, nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_cuda,
201       nbnxn_kernel_ElecRF_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecRF_VdwLJFsw_VF_cuda,
202       nbnxn_kernel_ElecRF_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_cuda,
203       nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_cuda },
204     { nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_VF_cuda,
205       nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_cuda,
206       nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_VF_cuda,
207       nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_VF_cuda },
208     { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_VF_cuda,
209       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_cuda,
210       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_VF_cuda,
211       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_VF_cuda },
212     { nbnxn_kernel_ElecEw_VdwLJ_VF_cuda, nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_cuda,
213       nbnxn_kernel_ElecEw_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecEw_VdwLJFsw_VF_cuda,
214       nbnxn_kernel_ElecEw_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_cuda,
215       nbnxn_kernel_ElecEw_VdwLJEwCombLB_VF_cuda },
216     { nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_cuda,
217       nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_cuda,
218       nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_cuda,
219       nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_VF_cuda }
220 };
221
222 /*! Force + pruning kernel function pointers. */
223 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_prune_ptr[c_numElecTypes][c_numVdwTypes] = {
224     { nbnxn_kernel_ElecCut_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecCut_VdwLJCombGeom_F_prune_cuda,
225       nbnxn_kernel_ElecCut_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecCut_VdwLJFsw_F_prune_cuda,
226       nbnxn_kernel_ElecCut_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombGeom_F_prune_cuda,
227       nbnxn_kernel_ElecCut_VdwLJEwCombLB_F_prune_cuda },
228     { nbnxn_kernel_ElecRF_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecRF_VdwLJCombGeom_F_prune_cuda,
229       nbnxn_kernel_ElecRF_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecRF_VdwLJFsw_F_prune_cuda,
230       nbnxn_kernel_ElecRF_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_prune_cuda,
231       nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_prune_cuda },
232     { nbnxn_kernel_ElecEwQSTab_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_F_prune_cuda,
233       nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_prune_cuda,
234       nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_F_prune_cuda,
235       nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_F_prune_cuda },
236     { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_prune_cuda,
237       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_F_prune_cuda,
238       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_F_prune_cuda,
239       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_prune_cuda,
240       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_F_prune_cuda,
241       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_F_prune_cuda },
242     { nbnxn_kernel_ElecEw_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecEw_VdwLJCombGeom_F_prune_cuda,
243       nbnxn_kernel_ElecEw_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecEw_VdwLJFsw_F_prune_cuda,
244       nbnxn_kernel_ElecEw_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_prune_cuda,
245       nbnxn_kernel_ElecEw_VdwLJEwCombLB_F_prune_cuda },
246     { nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_prune_cuda,
247       nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_prune_cuda,
248       nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_prune_cuda,
249       nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_F_prune_cuda }
250 };
251
252 /*! Force + energy + pruning kernel function pointers. */
253 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_prune_ptr[c_numElecTypes][c_numVdwTypes] = {
254     { nbnxn_kernel_ElecCut_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecCut_VdwLJCombGeom_VF_prune_cuda,
255       nbnxn_kernel_ElecCut_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecCut_VdwLJFsw_VF_prune_cuda,
256       nbnxn_kernel_ElecCut_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombGeom_VF_prune_cuda,
257       nbnxn_kernel_ElecCut_VdwLJEwCombLB_VF_prune_cuda },
258     { nbnxn_kernel_ElecRF_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_prune_cuda,
259       nbnxn_kernel_ElecRF_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecRF_VdwLJFsw_VF_prune_cuda,
260       nbnxn_kernel_ElecRF_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_prune_cuda,
261       nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_prune_cuda },
262     { nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_VF_prune_cuda,
263       nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_prune_cuda,
264       nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_VF_prune_cuda,
265       nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_VF_prune_cuda },
266     { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_prune_cuda,
267       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_VF_prune_cuda,
268       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_VF_prune_cuda,
269       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_prune_cuda,
270       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_prune_cuda,
271       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_VF_prune_cuda,
272       nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_VF_prune_cuda },
273     { nbnxn_kernel_ElecEw_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_prune_cuda,
274       nbnxn_kernel_ElecEw_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecEw_VdwLJFsw_VF_prune_cuda,
275       nbnxn_kernel_ElecEw_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_prune_cuda,
276       nbnxn_kernel_ElecEw_VdwLJEwCombLB_VF_prune_cuda },
277     { nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_prune_cuda,
278       nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_prune_cuda,
279       nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_prune_cuda,
280       nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_prune_cuda,
281       nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_VF_prune_cuda }
282 };
283
284 /*! Return a pointer to the kernel version to be executed at the current step. */
285 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(enum ElecType           elecType,
286                                                        enum VdwType            vdwType,
287                                                        bool                    bDoEne,
288                                                        bool                    bDoPrune,
289                                                        const DeviceInformation gmx_unused* deviceInfo)
290 {
291     const int elecTypeIdx = static_cast<int>(elecType);
292     const int vdwTypeIdx  = static_cast<int>(vdwType);
293
294     GMX_ASSERT(elecTypeIdx < c_numElecTypes,
295                "The electrostatics type requested is not implemented in the CUDA kernels.");
296     GMX_ASSERT(vdwTypeIdx < c_numVdwTypes,
297                "The VdW type requested is not implemented in the CUDA kernels.");
298
299     /* assert assumptions made by the kernels */
300     GMX_ASSERT(c_nbnxnGpuClusterSize * c_nbnxnGpuClusterSize / c_nbnxnGpuClusterpairSplit
301                        == deviceInfo->prop.warpSize,
302                "The CUDA kernels require the "
303                "cluster_size_i*cluster_size_j/nbnxn_gpu_clusterpair_split to match the warp size "
304                "of the architecture targeted.");
305
306     if (bDoEne)
307     {
308         if (bDoPrune)
309         {
310             return nb_kfunc_ener_prune_ptr[elecTypeIdx][vdwTypeIdx];
311         }
312         else
313         {
314             return nb_kfunc_ener_noprune_ptr[elecTypeIdx][vdwTypeIdx];
315         }
316     }
317     else
318     {
319         if (bDoPrune)
320         {
321             return nb_kfunc_noener_prune_ptr[elecTypeIdx][vdwTypeIdx];
322         }
323         else
324         {
325             return nb_kfunc_noener_noprune_ptr[elecTypeIdx][vdwTypeIdx];
326         }
327     }
328 }
329
330 /*! \brief Calculates the amount of shared memory required by the nonbonded kernel in use. */
331 static inline int calc_shmem_required_nonbonded(const int               num_threads_z,
332                                                 const DeviceInformation gmx_unused* deviceInfo,
333                                                 const NBParamGpu*                   nbp)
334 {
335     int shmem;
336
337     assert(deviceInfo);
338
339     /* size of shmem (force-buffers/xq/atom type preloading) */
340     /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
341     /* i-atom x+q in shared memory */
342     shmem = c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(float4);
343     /* cj in shared memory, for each warp separately */
344     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
345
346     if (nbp->vdwType == VdwType::CutCombGeom || nbp->vdwType == VdwType::CutCombLB)
347     {
348         /* i-atom LJ combination parameters in shared memory */
349         shmem += c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(float2);
350     }
351     else
352     {
353         /* i-atom types in shared memory */
354         shmem += c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(int);
355     }
356
357     return shmem;
358 }
359
360 /*! \brief Sync the nonlocal stream with dependent tasks in the local queue.
361  *
362  *  As the point where the local stream tasks can be considered complete happens
363  *  at the same call point where the nonlocal stream should be synced with the
364  *  the local, this function records the event if called with the local stream as
365  *  argument and inserts in the GPU stream a wait on the event on the nonlocal.
366  */
367 void nbnxnInsertNonlocalGpuDependency(const NbnxmGpu* nb, const InteractionLocality interactionLocality)
368 {
369     const DeviceStream& deviceStream = *nb->deviceStreams[interactionLocality];
370
371     /* When we get here all misc operations issued in the local stream as well as
372        the local xq H2D are done,
373        so we record that in the local stream and wait for it in the nonlocal one.
374        This wait needs to precede any PP tasks, bonded or nonbonded, that may
375        compute on interactions between local and nonlocal atoms.
376      */
377     if (nb->bUseTwoStreams)
378     {
379         if (interactionLocality == InteractionLocality::Local)
380         {
381             cudaError_t stat = cudaEventRecord(nb->misc_ops_and_local_H2D_done, deviceStream.stream());
382             CU_RET_ERR(stat, "cudaEventRecord on misc_ops_and_local_H2D_done failed");
383         }
384         else
385         {
386             cudaError_t stat =
387                     cudaStreamWaitEvent(deviceStream.stream(), nb->misc_ops_and_local_H2D_done, 0);
388             CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_and_local_H2D_done failed");
389         }
390     }
391 }
392
393 /*! \brief Launch asynchronously the xq buffer host to device copy. */
394 void gpu_copy_xq_to_gpu(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom, const AtomLocality atomLocality)
395 {
396     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
397
398     GMX_ASSERT(atomLocality == AtomLocality::Local || atomLocality == AtomLocality::NonLocal,
399                "Only local and non-local xq transfers are supported");
400
401     const InteractionLocality iloc = gpuAtomToInteractionLocality(atomLocality);
402
403     int adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
404
405     cu_atomdata_t*      adat         = nb->atdat;
406     gpu_plist*          plist        = nb->plist[iloc];
407     cu_timers_t*        t            = nb->timers;
408     const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
409
410     bool bDoTime = nb->bDoTime;
411
412     /* Don't launch the non-local H2D copy if there is no dependent
413        work to do: neither non-local nor other (e.g. bonded) work
414        to do that has as input the nbnxn coordaintes.
415        Doing the same for the local kernel is more complicated, since the
416        local part of the force array also depends on the non-local kernel.
417        So to avoid complicating the code and to reduce the risk of bugs,
418        we always call the local local x+q copy (and the rest of the local
419        work in nbnxn_gpu_launch_kernel().
420      */
421     if ((iloc == InteractionLocality::NonLocal) && !haveGpuShortRangeWork(*nb, iloc))
422     {
423         plist->haveFreshList = false;
424
425         return;
426     }
427
428     /* calculate the atom data index range based on locality */
429     if (atomLocality == AtomLocality::Local)
430     {
431         adat_begin = 0;
432         adat_len   = adat->natoms_local;
433     }
434     else
435     {
436         adat_begin = adat->natoms_local;
437         adat_len   = adat->natoms - adat->natoms_local;
438     }
439
440     /* HtoD x, q */
441     /* beginning of timed HtoD section */
442     if (bDoTime)
443     {
444         t->xf[atomLocality].nb_h2d.openTimingRegion(deviceStream);
445     }
446
447     static_assert(sizeof(adat->xq[0]) == sizeof(float4),
448                   "The size of the xyzq buffer element should be equal to the size of float4.");
449     copyToDeviceBuffer(&adat->xq, reinterpret_cast<const float4*>(nbatom->x().data()) + adat_begin,
450                        adat_begin, adat_len, deviceStream, GpuApiCallBehavior::Async, nullptr);
451
452     if (bDoTime)
453     {
454         t->xf[atomLocality].nb_h2d.closeTimingRegion(deviceStream);
455     }
456
457     /* When we get here all misc operations issued in the local stream as well as
458        the local xq H2D are done,
459        so we record that in the local stream and wait for it in the nonlocal one.
460        This wait needs to precede any PP tasks, bonded or nonbonded, that may
461        compute on interactions between local and nonlocal atoms.
462      */
463     nbnxnInsertNonlocalGpuDependency(nb, iloc);
464 }
465
466 /*! As we execute nonbonded workload in separate streams, before launching
467    the kernel we need to make sure that he following operations have completed:
468    - atomdata allocation and related H2D transfers (every nstlist step);
469    - pair list H2D transfer (every nstlist step);
470    - shift vector H2D transfer (every nstlist step);
471    - force (+shift force and energy) output clearing (every step).
472
473    These operations are issued in the local stream at the beginning of the step
474    and therefore always complete before the local kernel launch. The non-local
475    kernel is launched after the local on the same device/context hence it is
476    inherently scheduled after the operations in the local stream (including the
477    above "misc_ops") on pre-GK110 devices with single hardware queue, but on later
478    devices with multiple hardware queues the dependency needs to be enforced.
479    We use the misc_ops_and_local_H2D_done event to record the point where
480    the local x+q H2D (and all preceding) tasks are complete and synchronize
481    with this event in the non-local stream before launching the non-bonded kernel.
482  */
483 void gpu_launch_kernel(NbnxmGpu* nb, const gmx::StepWorkload& stepWork, const InteractionLocality iloc)
484 {
485     cu_atomdata_t*      adat         = nb->atdat;
486     NBParamGpu*         nbp          = nb->nbparam;
487     gpu_plist*          plist        = nb->plist[iloc];
488     cu_timers_t*        t            = nb->timers;
489     const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
490
491     bool bDoTime = nb->bDoTime;
492
493     /* Don't launch the non-local kernel if there is no work to do.
494        Doing the same for the local kernel is more complicated, since the
495        local part of the force array also depends on the non-local kernel.
496        So to avoid complicating the code and to reduce the risk of bugs,
497        we always call the local kernel, and later (not in
498        this function) the stream wait, local f copyback and the f buffer
499        clearing. All these operations, except for the local interaction kernel,
500        are needed for the non-local interactions. The skip of the local kernel
501        call is taken care of later in this function. */
502     if (canSkipNonbondedWork(*nb, iloc))
503     {
504         plist->haveFreshList = false;
505
506         return;
507     }
508
509     if (nbp->useDynamicPruning && plist->haveFreshList)
510     {
511         /* Prunes for rlistOuter and rlistInner, sets plist->haveFreshList=false
512            (TODO: ATM that's the way the timing accounting can distinguish between
513            separate prune kernel and combined force+prune, maybe we need a better way?).
514          */
515         gpu_launch_kernel_pruneonly(nb, iloc, 1);
516     }
517
518     if (plist->nsci == 0)
519     {
520         /* Don't launch an empty local kernel (not allowed with CUDA) */
521         return;
522     }
523
524     /* beginning of timed nonbonded calculation section */
525     if (bDoTime)
526     {
527         t->interaction[iloc].nb_k.openTimingRegion(deviceStream);
528     }
529
530     /* Kernel launch config:
531      * - The thread block dimensions match the size of i-clusters, j-clusters,
532      *   and j-cluster concurrency, in x, y, and z, respectively.
533      * - The 1D block-grid contains as many blocks as super-clusters.
534      */
535     int num_threads_z = 1;
536     if (nb->deviceContext_->deviceInfo().prop.major == 3 && nb->deviceContext_->deviceInfo().prop.minor == 7)
537     {
538         num_threads_z = 2;
539     }
540     int nblock = calc_nb_kernel_nblock(plist->nsci, &nb->deviceContext_->deviceInfo());
541
542
543     KernelLaunchConfig config;
544     config.blockSize[0] = c_clSize;
545     config.blockSize[1] = c_clSize;
546     config.blockSize[2] = num_threads_z;
547     config.gridSize[0]  = nblock;
548     config.sharedMemorySize =
549             calc_shmem_required_nonbonded(num_threads_z, &nb->deviceContext_->deviceInfo(), nbp);
550
551     if (debug)
552     {
553         fprintf(debug,
554                 "Non-bonded GPU launch configuration:\n\tThread block: %zux%zux%zu\n\t"
555                 "\tGrid: %zux%zu\n\t#Super-clusters/clusters: %d/%d (%d)\n"
556                 "\tShMem: %zu\n",
557                 config.blockSize[0], config.blockSize[1], config.blockSize[2], config.gridSize[0],
558                 config.gridSize[1], plist->nsci * c_nbnxnGpuNumClusterPerSupercluster,
559                 c_nbnxnGpuNumClusterPerSupercluster, plist->na_c, config.sharedMemorySize);
560     }
561
562     auto*      timingEvent = bDoTime ? t->interaction[iloc].nb_k.fetchNextEvent() : nullptr;
563     const auto kernel =
564             select_nbnxn_kernel(nbp->elecType, nbp->vdwType, stepWork.computeEnergy,
565                                 (plist->haveFreshList && !nb->timers->interaction[iloc].didPrune),
566                                 &nb->deviceContext_->deviceInfo());
567     const auto kernelArgs =
568             prepareGpuKernelArguments(kernel, config, adat, nbp, plist, &stepWork.computeVirial);
569     launchGpuKernel(kernel, config, deviceStream, timingEvent, "k_calc_nb", kernelArgs);
570
571     if (bDoTime)
572     {
573         t->interaction[iloc].nb_k.closeTimingRegion(deviceStream);
574     }
575
576     if (GMX_NATIVE_WINDOWS)
577     {
578         /* Windows: force flushing WDDM queue */
579         cudaStreamQuery(deviceStream.stream());
580     }
581 }
582
583 /*! Calculates the amount of shared memory required by the CUDA kernel in use. */
584 static inline int calc_shmem_required_prune(const int num_threads_z)
585 {
586     int shmem;
587
588     /* i-atom x in shared memory */
589     shmem = c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(float4);
590     /* cj in shared memory, for each warp separately */
591     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
592
593     return shmem;
594 }
595
596 void gpu_launch_kernel_pruneonly(NbnxmGpu* nb, const InteractionLocality iloc, const int numParts)
597 {
598     cu_atomdata_t*      adat         = nb->atdat;
599     NBParamGpu*         nbp          = nb->nbparam;
600     gpu_plist*          plist        = nb->plist[iloc];
601     cu_timers_t*        t            = nb->timers;
602     const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
603
604     bool bDoTime = nb->bDoTime;
605
606     if (plist->haveFreshList)
607     {
608         GMX_ASSERT(numParts == 1, "With first pruning we expect 1 part");
609
610         /* Set rollingPruningNumParts to signal that it is not set */
611         plist->rollingPruningNumParts = 0;
612         plist->rollingPruningPart     = 0;
613     }
614     else
615     {
616         if (plist->rollingPruningNumParts == 0)
617         {
618             plist->rollingPruningNumParts = numParts;
619         }
620         else
621         {
622             GMX_ASSERT(numParts == plist->rollingPruningNumParts,
623                        "It is not allowed to change numParts in between list generation steps");
624         }
625     }
626
627     /* Use a local variable for part and update in plist, so we can return here
628      * without duplicating the part increment code.
629      */
630     int part = plist->rollingPruningPart;
631
632     plist->rollingPruningPart++;
633     if (plist->rollingPruningPart >= plist->rollingPruningNumParts)
634     {
635         plist->rollingPruningPart = 0;
636     }
637
638     /* Compute the number of list entries to prune in this pass */
639     int numSciInPart = (plist->nsci - part) / numParts;
640
641     /* Don't launch the kernel if there is no work to do (not allowed with CUDA) */
642     if (numSciInPart <= 0)
643     {
644         plist->haveFreshList = false;
645
646         return;
647     }
648
649     GpuRegionTimer* timer = nullptr;
650     if (bDoTime)
651     {
652         timer = &(plist->haveFreshList ? t->interaction[iloc].prune_k : t->interaction[iloc].rollingPrune_k);
653     }
654
655     /* beginning of timed prune calculation section */
656     if (bDoTime)
657     {
658         timer->openTimingRegion(deviceStream);
659     }
660
661     /* Kernel launch config:
662      * - The thread block dimensions match the size of i-clusters, j-clusters,
663      *   and j-cluster concurrency, in x, y, and z, respectively.
664      * - The 1D block-grid contains as many blocks as super-clusters.
665      */
666     int num_threads_z = c_cudaPruneKernelJ4Concurrency;
667     int nblock        = calc_nb_kernel_nblock(numSciInPart, &nb->deviceContext_->deviceInfo());
668     KernelLaunchConfig config;
669     config.blockSize[0]     = c_clSize;
670     config.blockSize[1]     = c_clSize;
671     config.blockSize[2]     = num_threads_z;
672     config.gridSize[0]      = nblock;
673     config.sharedMemorySize = calc_shmem_required_prune(num_threads_z);
674
675     if (debug)
676     {
677         fprintf(debug,
678                 "Pruning GPU kernel launch configuration:\n\tThread block: %zux%zux%zu\n\t"
679                 "\tGrid: %zux%zu\n\t#Super-clusters/clusters: %d/%d (%d)\n"
680                 "\tShMem: %zu\n",
681                 config.blockSize[0], config.blockSize[1], config.blockSize[2], config.gridSize[0],
682                 config.gridSize[1], numSciInPart * c_nbnxnGpuNumClusterPerSupercluster,
683                 c_nbnxnGpuNumClusterPerSupercluster, plist->na_c, config.sharedMemorySize);
684     }
685
686     auto*          timingEvent  = bDoTime ? timer->fetchNextEvent() : nullptr;
687     constexpr char kernelName[] = "k_pruneonly";
688     const auto     kernel =
689             plist->haveFreshList ? nbnxn_kernel_prune_cuda<true> : nbnxn_kernel_prune_cuda<false>;
690     const auto kernelArgs = prepareGpuKernelArguments(kernel, config, adat, nbp, plist, &numParts, &part);
691     launchGpuKernel(kernel, config, deviceStream, timingEvent, kernelName, kernelArgs);
692
693     /* TODO: consider a more elegant way to track which kernel has been called
694        (combined or separate 1st pass prune, rolling prune). */
695     if (plist->haveFreshList)
696     {
697         plist->haveFreshList = false;
698         /* Mark that pruning has been done */
699         nb->timers->interaction[iloc].didPrune = true;
700     }
701     else
702     {
703         /* Mark that rolling pruning has been done */
704         nb->timers->interaction[iloc].didRollingPrune = true;
705     }
706
707     if (bDoTime)
708     {
709         timer->closeTimingRegion(deviceStream);
710     }
711
712     if (GMX_NATIVE_WINDOWS)
713     {
714         /* Windows: force flushing WDDM queue */
715         cudaStreamQuery(deviceStream.stream());
716     }
717 }
718
719 void gpu_launch_cpyback(NbnxmGpu*                nb,
720                         nbnxn_atomdata_t*        nbatom,
721                         const gmx::StepWorkload& stepWork,
722                         const AtomLocality       atomLocality)
723 {
724     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
725
726     cudaError_t stat;
727     int         adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
728
729     /* determine interaction locality from atom locality */
730     const InteractionLocality iloc = gpuAtomToInteractionLocality(atomLocality);
731
732     /* extract the data */
733     cu_atomdata_t*      adat         = nb->atdat;
734     cu_timers_t*        t            = nb->timers;
735     bool                bDoTime      = nb->bDoTime;
736     const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
737
738     /* don't launch non-local copy-back if there was no non-local work to do */
739     if ((iloc == InteractionLocality::NonLocal) && !haveGpuShortRangeWork(*nb, iloc))
740     {
741         return;
742     }
743
744     getGpuAtomRange(adat, atomLocality, &adat_begin, &adat_len);
745
746     /* beginning of timed D2H section */
747     if (bDoTime)
748     {
749         t->xf[atomLocality].nb_d2h.openTimingRegion(deviceStream);
750     }
751
752     /* With DD the local D2H transfer can only start after the non-local
753        kernel has finished. */
754     if (iloc == InteractionLocality::Local && nb->bUseTwoStreams)
755     {
756         stat = cudaStreamWaitEvent(deviceStream.stream(), nb->nonlocal_done, 0);
757         CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
758     }
759
760     /* DtoH f
761      * Skip if buffer ops / reduction is offloaded to the GPU.
762      */
763     if (!stepWork.useGpuFBufferOps)
764     {
765         static_assert(
766                 sizeof(adat->f[0]) == sizeof(float3),
767                 "The size of the force buffer element should be equal to the size of float3.");
768         copyFromDeviceBuffer(reinterpret_cast<float3*>(nbatom->out[0].f.data()) + adat_begin, &adat->f,
769                              adat_begin, adat_len, deviceStream, GpuApiCallBehavior::Async, nullptr);
770     }
771
772     /* After the non-local D2H is launched the nonlocal_done event can be
773        recorded which signals that the local D2H can proceed. This event is not
774        placed after the non-local kernel because we want the non-local data
775        back first. */
776     if (iloc == InteractionLocality::NonLocal)
777     {
778         stat = cudaEventRecord(nb->nonlocal_done, deviceStream.stream());
779         CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
780     }
781
782     /* only transfer energies in the local stream */
783     if (iloc == InteractionLocality::Local)
784     {
785         /* DtoH fshift when virial is needed */
786         if (stepWork.computeVirial)
787         {
788             static_assert(sizeof(nb->nbst.fshift[0]) == sizeof(adat->fshift[0]),
789                           "Sizes of host- and device-side shift vectors should be the same.");
790             copyFromDeviceBuffer(nb->nbst.fshift, &adat->fshift, 0, SHIFTS, deviceStream,
791                                  GpuApiCallBehavior::Async, nullptr);
792         }
793
794         /* DtoH energies */
795         if (stepWork.computeEnergy)
796         {
797             static_assert(sizeof(nb->nbst.e_lj[0]) == sizeof(adat->e_lj[0]),
798                           "Sizes of host- and device-side LJ energy terms should be the same.");
799             copyFromDeviceBuffer(nb->nbst.e_lj, &adat->e_lj, 0, 1, deviceStream,
800                                  GpuApiCallBehavior::Async, nullptr);
801             static_assert(sizeof(nb->nbst.e_el[0]) == sizeof(adat->e_el[0]),
802                           "Sizes of host- and device-side electrostatic energy terms should be the "
803                           "same.");
804             copyFromDeviceBuffer(nb->nbst.e_el, &adat->e_el, 0, 1, deviceStream,
805                                  GpuApiCallBehavior::Async, nullptr);
806         }
807     }
808
809     if (bDoTime)
810     {
811         t->xf[atomLocality].nb_d2h.closeTimingRegion(deviceStream);
812     }
813 }
814
815 void cuda_set_cacheconfig()
816 {
817     cudaError_t stat;
818
819     for (int i = 0; i < c_numElecTypes; i++)
820     {
821         for (int j = 0; j < c_numVdwTypes; j++)
822         {
823             /* Default kernel 32/32 kB Shared/L1 */
824             cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferEqual);
825             cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
826             cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferEqual);
827             stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
828             CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");
829         }
830     }
831 }
832
833 /* X buffer operations on GPU: performs conversion from rvec to nb format. */
834 void nbnxn_gpu_x_to_nbat_x(const Nbnxm::Grid&        grid,
835                            bool                      setFillerCoords,
836                            NbnxmGpu*                 nb,
837                            DeviceBuffer<gmx::RVec>   d_x,
838                            GpuEventSynchronizer*     xReadyOnDevice,
839                            const Nbnxm::AtomLocality locality,
840                            int                       gridId,
841                            int                       numColumnsMax)
842 {
843     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
844
845     cu_atomdata_t* adat = nb->atdat;
846
847     const int                  numColumns      = grid.numColumns();
848     const int                  cellOffset      = grid.cellOffset();
849     const int                  numAtomsPerCell = grid.numAtomsPerCell();
850     Nbnxm::InteractionLocality interactionLoc  = gpuAtomToInteractionLocality(locality);
851
852     const DeviceStream& deviceStream = *nb->deviceStreams[interactionLoc];
853
854     int numAtoms = grid.srcAtomEnd() - grid.srcAtomBegin();
855     // avoid empty kernel launch, skip to inserting stream dependency
856     if (numAtoms != 0)
857     {
858         // TODO: This will only work with CUDA
859         GMX_ASSERT(d_x, "Need a valid device pointer");
860
861         // ensure that coordinates are ready on the device before launching the kernel
862         GMX_ASSERT(xReadyOnDevice, "Need a valid GpuEventSynchronizer object");
863         xReadyOnDevice->enqueueWaitEvent(deviceStream);
864
865         KernelLaunchConfig config;
866         config.blockSize[0] = c_bufOpsThreadsPerBlock;
867         config.blockSize[1] = 1;
868         config.blockSize[2] = 1;
869         config.gridSize[0] = (grid.numCellsColumnMax() * numAtomsPerCell + c_bufOpsThreadsPerBlock - 1)
870                              / c_bufOpsThreadsPerBlock;
871         config.gridSize[1] = numColumns;
872         config.gridSize[2] = 1;
873         GMX_ASSERT(config.gridSize[0] > 0,
874                    "Can not have empty grid, early return above avoids this");
875         config.sharedMemorySize = 0;
876
877         auto kernelFn = setFillerCoords ? nbnxn_gpu_x_to_nbat_x_kernel<true>
878                                         : nbnxn_gpu_x_to_nbat_x_kernel<false>;
879         float4*    d_xq          = adat->xq;
880         float3*    d_xFloat3     = asFloat3(d_x);
881         const int* d_atomIndices = nb->atomIndices;
882         const int* d_cxy_na      = &nb->cxy_na[numColumnsMax * gridId];
883         const int* d_cxy_ind     = &nb->cxy_ind[numColumnsMax * gridId];
884         const auto kernelArgs    = prepareGpuKernelArguments(kernelFn, config, &numColumns, &d_xq,
885                                                           &d_xFloat3, &d_atomIndices, &d_cxy_na,
886                                                           &d_cxy_ind, &cellOffset, &numAtomsPerCell);
887         launchGpuKernel(kernelFn, config, deviceStream, nullptr, "XbufferOps", kernelArgs);
888     }
889
890     // TODO: note that this is not necessary when there astreamre no local atoms, that is:
891     // (numAtoms == 0 && interactionLoc == InteractionLocality::Local)
892     // but for now we avoid that optimization
893     nbnxnInsertNonlocalGpuDependency(nb, interactionLoc);
894 }
895
896 void* getGpuForces(NbnxmGpu* nb)
897 {
898     return nb->atdat->f;
899 }
900
901 } // namespace Nbnxm