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