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