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