05c0278abc9c71878138f137e0efc2d4819d07fa
[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/gpu_utils/gpueventsynchronizer.cuh"
58 #include "gromacs/gpu_utils/vectype_ops.cuh"
59 #include "gromacs/mdlib/ppforceworkload.h"
60 #include "gromacs/nbnxm/atomdata.h"
61 #include "gromacs/nbnxm/gpu_common.h"
62 #include "gromacs/nbnxm/gpu_common_utils.h"
63 #include "gromacs/nbnxm/gpu_data_mgmt.h"
64 #include "gromacs/nbnxm/grid.h"
65 #include "gromacs/nbnxm/nbnxm.h"
66 #include "gromacs/nbnxm/pairlist.h"
67 #include "gromacs/nbnxm/cuda/nbnxm_buffer_ops_kernels.cuh"
68 #include "gromacs/timing/gpu_timing.h"
69 #include "gromacs/utility/cstringutil.h"
70 #include "gromacs/utility/gmxassert.h"
71
72 #include "nbnxm_cuda_types.h"
73
74 /***** The kernel declarations/definitions come here *****/
75
76 /* Top-level kernel declaration generation: will generate through multiple
77  * inclusion the following flavors for all kernel declarations:
78  * - force-only output;
79  * - force and energy output;
80  * - force-only with pair list pruning;
81  * - force and energy output with pair list pruning.
82  */
83 #define FUNCTION_DECLARATION_ONLY
84 /** Force only **/
85 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernels.cuh"
86 /** Force & energy **/
87 #define CALC_ENERGIES
88 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernels.cuh"
89 #undef CALC_ENERGIES
90
91 /*** Pair-list pruning kernels ***/
92 /** Force only **/
93 #define PRUNE_NBL
94 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernels.cuh"
95 /** Force & energy **/
96 #define CALC_ENERGIES
97 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernels.cuh"
98 #undef CALC_ENERGIES
99 #undef PRUNE_NBL
100
101 /* Prune-only kernels */
102 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernel_pruneonly.cuh"
103 #undef FUNCTION_DECLARATION_ONLY
104
105 /* Now generate the function definitions if we are using a single compilation unit. */
106 #if GMX_CUDA_NB_SINGLE_COMPILATION_UNIT
107 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernel_F_noprune.cu"
108 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernel_F_prune.cu"
109 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernel_VF_noprune.cu"
110 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernel_VF_prune.cu"
111 #include "gromacs/nbnxm/cuda/nbnxm_cuda_kernel_pruneonly.cu"
112 #endif /* GMX_CUDA_NB_SINGLE_COMPILATION_UNIT */
113
114 namespace Nbnxm
115 {
116
117 //! Number of CUDA threads in a block
118 //TODO Optimize this through experimentation
119 constexpr static int c_bufOpsThreadsPerBlock = 128;
120
121 /*! Nonbonded kernel function pointer type */
122 typedef void (*nbnxn_cu_kfunc_ptr_t)(const cu_atomdata_t,
123                                      const cu_nbparam_t,
124                                      const cu_plist_t,
125                                      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 gmx_device_info_t *dinfo)
131 {
132     int max_grid_x_size;
133
134     assert(dinfo);
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 = dinfo->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, "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)!", nwork_units, max_grid_x_size);
147     }
148
149     return nwork_units;
150 }
151
152
153 /* Constant arrays listing all kernel function pointers and enabling selection
154    of a kernel in an elegant manner. */
155
156 /*! Pointers to the non-bonded kernels organized in 2-dim arrays by:
157  *  electrostatics and VDW type.
158  *
159  *  Note that the row- and column-order of function pointers has to match the
160  *  order of corresponding enumerated electrostatics and vdw types, resp.,
161  *  defined in nbnxn_cuda_types.h.
162  */
163
164 /*! Force-only kernel function pointers. */
165 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_noprune_ptr[eelCuNR][evdwCuNR] =
166 {
167     { 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            },
168     { 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             },
169     { 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        },
170     { 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 },
171     { 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             },
172     { 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      }
173 };
174
175 /*! Force + energy kernel function pointers. */
176 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_noprune_ptr[eelCuNR][evdwCuNR] =
177 {
178     { 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            },
179     { 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             },
180     { 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        },
181     { 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 },
182     { 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             },
183     { 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      }
184 };
185
186 /*! Force + pruning kernel function pointers. */
187 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_prune_ptr[eelCuNR][evdwCuNR] =
188 {
189     { 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             },
190     { 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              },
191     { 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         },
192     { 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  },
193     { 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              },
194     { 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       }
195 };
196
197 /*! Force + energy + pruning kernel function pointers. */
198 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_prune_ptr[eelCuNR][evdwCuNR] =
199 {
200     { 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            },
201     { 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             },
202     { 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        },
203     { 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 },
204     { 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             },
205     { 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      }
206 };
207
208 /*! Return a pointer to the kernel version to be executed at the current step. */
209 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(int                                  eeltype,
210                                                        int                                  evdwtype,
211                                                        bool                                 bDoEne,
212                                                        bool                                 bDoPrune,
213                                                        const gmx_device_info_t gmx_unused  *devInfo)
214 {
215     nbnxn_cu_kfunc_ptr_t res;
216
217     GMX_ASSERT(eeltype < eelCuNR,
218                "The electrostatics type requested is not implemented in the CUDA kernels.");
219     GMX_ASSERT(evdwtype < evdwCuNR,
220                "The VdW type requested is not implemented in the CUDA kernels.");
221
222     /* assert assumptions made by the kernels */
223     GMX_ASSERT(c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit == devInfo->prop.warpSize,
224                "The CUDA kernels require the cluster_size_i*cluster_size_j/nbnxn_gpu_clusterpair_split to match the warp size of the architecture targeted.");
225
226     if (bDoEne)
227     {
228         if (bDoPrune)
229         {
230             res = nb_kfunc_ener_prune_ptr[eeltype][evdwtype];
231         }
232         else
233         {
234             res = nb_kfunc_ener_noprune_ptr[eeltype][evdwtype];
235         }
236     }
237     else
238     {
239         if (bDoPrune)
240         {
241             res = nb_kfunc_noener_prune_ptr[eeltype][evdwtype];
242         }
243         else
244         {
245             res = nb_kfunc_noener_noprune_ptr[eeltype][evdwtype];
246         }
247     }
248
249     return res;
250 }
251
252 /*! \brief Calculates the amount of shared memory required by the nonbonded kernel in use. */
253 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)
254 {
255     int shmem;
256
257     assert(dinfo);
258
259     /* size of shmem (force-buffers/xq/atom type preloading) */
260     /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
261     /* i-atom x+q in shared memory */
262     shmem  = c_numClPerSupercl * c_clSize * sizeof(float4);
263     /* cj in shared memory, for each warp separately */
264     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
265
266     if (nbp->vdwtype == evdwCuCUTCOMBGEOM ||
267         nbp->vdwtype == evdwCuCUTCOMBLB)
268     {
269         /* i-atom LJ combination parameters in shared memory */
270         shmem += c_numClPerSupercl * c_clSize * sizeof(float2);
271     }
272     else
273     {
274         /* i-atom types in shared memory */
275         shmem += c_numClPerSupercl * c_clSize * sizeof(int);
276     }
277
278     return shmem;
279 }
280
281 /*! \brief Sync the nonlocal stream with dependent tasks in the local queue.
282  *
283  *  As the point where the local stream tasks can be considered complete happens
284  *  at the same call point where the nonlocal stream should be synced with the
285  *  the local, this function records the event if called with the local stream as
286  *  argument and inserts in the GPU stream a wait on the event on the nonlocal.
287  */
288 void nbnxnInsertNonlocalGpuDependency(const gmx_nbnxn_cuda_t   *nb,
289                                       const InteractionLocality interactionLocality)
290 {
291     cudaStream_t stream  = nb->stream[interactionLocality];
292
293     /* When we get here all misc operations issued in the local stream as well as
294        the local xq H2D are done,
295        so we record that in the local stream and wait for it in the nonlocal one.
296        This wait needs to precede any PP tasks, bonded or nonbonded, that may
297        compute on interactions between local and nonlocal atoms.
298      */
299     if (nb->bUseTwoStreams)
300     {
301         if (interactionLocality == InteractionLocality::Local)
302         {
303             cudaError_t stat = cudaEventRecord(nb->misc_ops_and_local_H2D_done, stream);
304             CU_RET_ERR(stat, "cudaEventRecord on misc_ops_and_local_H2D_done failed");
305         }
306         else
307         {
308             cudaError_t stat = cudaStreamWaitEvent(stream, nb->misc_ops_and_local_H2D_done, 0);
309             CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_and_local_H2D_done failed");
310         }
311     }
312 }
313
314 /*! \brief Launch asynchronously the xq buffer host to device copy. */
315 void gpu_copy_xq_to_gpu(gmx_nbnxn_cuda_t       *nb,
316                         const nbnxn_atomdata_t *nbatom,
317                         const AtomLocality      atomLocality)
318 {
319     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
320
321     GMX_ASSERT(atomLocality == AtomLocality::Local || atomLocality == AtomLocality::NonLocal,
322                "Only local and non-local xq transfers are supported");
323
324     const InteractionLocality iloc = gpuAtomToInteractionLocality(atomLocality);
325
326     int                       adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
327
328     cu_atomdata_t            *adat    = nb->atdat;
329     cu_plist_t               *plist   = nb->plist[iloc];
330     cu_timers_t              *t       = nb->timers;
331     cudaStream_t              stream  = nb->stream[iloc];
332
333     bool                      bDoTime     = nb->bDoTime;
334
335     /* Don't launch the non-local H2D copy if there is no dependent
336        work to do: neither non-local nor other (e.g. bonded) work
337        to do that has as input the nbnxn coordaintes.
338        Doing the same for the local kernel is more complicated, since the
339        local part of the force array also depends on the non-local kernel.
340        So to avoid complicating the code and to reduce the risk of bugs,
341        we always call the local local x+q copy (and the rest of the local
342        work in nbnxn_gpu_launch_kernel().
343      */
344     if ((iloc == InteractionLocality::NonLocal) && !haveGpuShortRangeWork(*nb, iloc))
345     {
346         plist->haveFreshList = false;
347
348         return;
349     }
350
351     /* calculate the atom data index range based on locality */
352     if (atomLocality == AtomLocality::Local)
353     {
354         adat_begin  = 0;
355         adat_len    = adat->natoms_local;
356     }
357     else
358     {
359         adat_begin  = adat->natoms_local;
360         adat_len    = adat->natoms - adat->natoms_local;
361     }
362
363     /* HtoD x, q */
364     /* beginning of timed HtoD section */
365     if (bDoTime)
366     {
367         t->xf[atomLocality].nb_h2d.openTimingRegion(stream);
368     }
369
370     cu_copy_H2D_async(adat->xq + adat_begin, static_cast<const void *>(nbatom->x().data() + adat_begin * 4),
371                       adat_len * sizeof(*adat->xq), stream);
372
373     if (bDoTime)
374     {
375         t->xf[atomLocality].nb_h2d.closeTimingRegion(stream);
376     }
377
378     /* When we get here all misc operations issued in the local stream as well as
379        the local xq H2D are done,
380        so we record that in the local stream and wait for it in the nonlocal one.
381        This wait needs to precede any PP tasks, bonded or nonbonded, that may
382        compute on interactions between local and nonlocal atoms.
383      */
384     nbnxnInsertNonlocalGpuDependency(nb, iloc);
385 }
386
387 /*! As we execute nonbonded workload in separate streams, before launching
388    the kernel we need to make sure that he following operations have completed:
389    - atomdata allocation and related H2D transfers (every nstlist step);
390    - pair list H2D transfer (every nstlist step);
391    - shift vector H2D transfer (every nstlist step);
392    - force (+shift force and energy) output clearing (every step).
393
394    These operations are issued in the local stream at the beginning of the step
395    and therefore always complete before the local kernel launch. The non-local
396    kernel is launched after the local on the same device/context hence it is
397    inherently scheduled after the operations in the local stream (including the
398    above "misc_ops") on pre-GK110 devices with single hardware queue, but on later
399    devices with multiple hardware queues the dependency needs to be enforced.
400    We use the misc_ops_and_local_H2D_done event to record the point where
401    the local x+q H2D (and all preceding) tasks are complete and synchronize
402    with this event in the non-local stream before launching the non-bonded kernel.
403  */
404 void gpu_launch_kernel(gmx_nbnxn_cuda_t          *nb,
405                        const gmx::ForceFlags     &forceFlags,
406                        const InteractionLocality  iloc)
407 {
408     cu_atomdata_t       *adat    = nb->atdat;
409     cu_nbparam_t        *nbp     = nb->nbparam;
410     cu_plist_t          *plist   = nb->plist[iloc];
411     cu_timers_t         *t       = nb->timers;
412     cudaStream_t         stream  = nb->stream[iloc];
413
414     bool                 bDoTime     = nb->bDoTime;
415
416     /* Don't launch the non-local kernel if there is no work to do.
417        Doing the same for the local kernel is more complicated, since the
418        local part of the force array also depends on the non-local kernel.
419        So to avoid complicating the code and to reduce the risk of bugs,
420        we always call the local kernel, and later (not in
421        this function) the stream wait, local f copyback and the f buffer
422        clearing. All these operations, except for the local interaction kernel,
423        are needed for the non-local interactions. The skip of the local kernel
424        call is taken care of later in this function. */
425     if (canSkipNonbondedWork(*nb, iloc))
426     {
427         plist->haveFreshList = false;
428
429         return;
430     }
431
432     if (nbp->useDynamicPruning && plist->haveFreshList)
433     {
434         /* Prunes for rlistOuter and rlistInner, sets plist->haveFreshList=false
435            (TODO: ATM that's the way the timing accounting can distinguish between
436            separate prune kernel and combined force+prune, maybe we need a better way?).
437          */
438         gpu_launch_kernel_pruneonly(nb, iloc, 1);
439     }
440
441     if (plist->nsci == 0)
442     {
443         /* Don't launch an empty local kernel (not allowed with CUDA) */
444         return;
445     }
446
447     /* beginning of timed nonbonded calculation section */
448     if (bDoTime)
449     {
450         t->interaction[iloc].nb_k.openTimingRegion(stream);
451     }
452
453     /* Kernel launch config:
454      * - The thread block dimensions match the size of i-clusters, j-clusters,
455      *   and j-cluster concurrency, in x, y, and z, respectively.
456      * - The 1D block-grid contains as many blocks as super-clusters.
457      */
458     int num_threads_z = 1;
459     if (nb->dev_info->prop.major == 3 && nb->dev_info->prop.minor == 7)
460     {
461         num_threads_z = 2;
462     }
463     int nblock    = calc_nb_kernel_nblock(plist->nsci, nb->dev_info);
464
465
466     KernelLaunchConfig config;
467     config.blockSize[0]     = c_clSize;
468     config.blockSize[1]     = c_clSize;
469     config.blockSize[2]     = num_threads_z;
470     config.gridSize[0]      = nblock;
471     config.sharedMemorySize = calc_shmem_required_nonbonded(num_threads_z, nb->dev_info, nbp);
472     config.stream           = stream;
473
474     if (debug)
475     {
476         fprintf(debug, "Non-bonded GPU launch configuration:\n\tThread block: %zux%zux%zu\n\t"
477                 "\tGrid: %zux%zu\n\t#Super-clusters/clusters: %d/%d (%d)\n"
478                 "\tShMem: %zu\n",
479                 config.blockSize[0], config.blockSize[1], config.blockSize[2],
480                 config.gridSize[0], config.gridSize[1], plist->nsci*c_numClPerSupercl,
481                 c_numClPerSupercl, plist->na_c,
482                 config.sharedMemorySize);
483     }
484
485     auto       *timingEvent = bDoTime ? t->interaction[iloc].nb_k.fetchNextEvent() : nullptr;
486     const auto  kernel      = select_nbnxn_kernel(nbp->eeltype,
487                                                   nbp->vdwtype,
488                                                   forceFlags.computeEnergy,
489                                                   (plist->haveFreshList && !nb->timers->interaction[iloc].didPrune),
490                                                   nb->dev_info);
491     const auto kernelArgs  = prepareGpuKernelArguments(kernel, config, adat, nbp, plist, &forceFlags.computeVirial);
492     launchGpuKernel(kernel, config, timingEvent, "k_calc_nb", kernelArgs);
493
494     if (bDoTime)
495     {
496         t->interaction[iloc].nb_k.closeTimingRegion(stream);
497     }
498
499     if (GMX_NATIVE_WINDOWS)
500     {
501         /* Windows: force flushing WDDM queue */
502         cudaStreamQuery(stream);
503     }
504 }
505
506 /*! Calculates the amount of shared memory required by the CUDA kernel in use. */
507 static inline int calc_shmem_required_prune(const int num_threads_z)
508 {
509     int shmem;
510
511     /* i-atom x in shared memory */
512     shmem  = c_numClPerSupercl * c_clSize * sizeof(float4);
513     /* cj in shared memory, for each warp separately */
514     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
515
516     return shmem;
517 }
518
519 void gpu_launch_kernel_pruneonly(gmx_nbnxn_cuda_t          *nb,
520                                  const InteractionLocality  iloc,
521                                  const int                  numParts)
522 {
523     cu_atomdata_t       *adat    = nb->atdat;
524     cu_nbparam_t        *nbp     = nb->nbparam;
525     cu_plist_t          *plist   = nb->plist[iloc];
526     cu_timers_t         *t       = nb->timers;
527     cudaStream_t         stream  = nb->stream[iloc];
528
529     bool                 bDoTime = nb->bDoTime;
530
531     if (plist->haveFreshList)
532     {
533         GMX_ASSERT(numParts == 1, "With first pruning we expect 1 part");
534
535         /* Set rollingPruningNumParts to signal that it is not set */
536         plist->rollingPruningNumParts = 0;
537         plist->rollingPruningPart     = 0;
538     }
539     else
540     {
541         if (plist->rollingPruningNumParts == 0)
542         {
543             plist->rollingPruningNumParts = numParts;
544         }
545         else
546         {
547             GMX_ASSERT(numParts == plist->rollingPruningNumParts, "It is not allowed to change numParts in between list generation steps");
548         }
549     }
550
551     /* Use a local variable for part and update in plist, so we can return here
552      * without duplicating the part increment code.
553      */
554     int part = plist->rollingPruningPart;
555
556     plist->rollingPruningPart++;
557     if (plist->rollingPruningPart >= plist->rollingPruningNumParts)
558     {
559         plist->rollingPruningPart = 0;
560     }
561
562     /* Compute the number of list entries to prune in this pass */
563     int numSciInPart = (plist->nsci - part)/numParts;
564
565     /* Don't launch the kernel if there is no work to do (not allowed with CUDA) */
566     if (numSciInPart <= 0)
567     {
568         plist->haveFreshList = false;
569
570         return;
571     }
572
573     GpuRegionTimer *timer = nullptr;
574     if (bDoTime)
575     {
576         timer = &(plist->haveFreshList ? t->interaction[iloc].prune_k : t->interaction[iloc].rollingPrune_k);
577     }
578
579     /* beginning of timed prune calculation section */
580     if (bDoTime)
581     {
582         timer->openTimingRegion(stream);
583     }
584
585     /* Kernel launch config:
586      * - The thread block dimensions match the size of i-clusters, j-clusters,
587      *   and j-cluster concurrency, in x, y, and z, respectively.
588      * - The 1D block-grid contains as many blocks as super-clusters.
589      */
590     int                num_threads_z  = c_cudaPruneKernelJ4Concurrency;
591     int                nblock         = calc_nb_kernel_nblock(numSciInPart, nb->dev_info);
592     KernelLaunchConfig config;
593     config.blockSize[0]     = c_clSize;
594     config.blockSize[1]     = c_clSize;
595     config.blockSize[2]     = num_threads_z;
596     config.gridSize[0]      = nblock;
597     config.sharedMemorySize = calc_shmem_required_prune(num_threads_z);
598     config.stream           = stream;
599
600     if (debug)
601     {
602         fprintf(debug, "Pruning GPU kernel launch configuration:\n\tThread block: %zux%zux%zu\n\t"
603                 "\tGrid: %zux%zu\n\t#Super-clusters/clusters: %d/%d (%d)\n"
604                 "\tShMem: %zu\n",
605                 config.blockSize[0], config.blockSize[1], config.blockSize[2],
606                 config.gridSize[0], config.gridSize[1], numSciInPart*c_numClPerSupercl,
607                 c_numClPerSupercl, plist->na_c,
608                 config.sharedMemorySize);
609     }
610
611     auto          *timingEvent  = bDoTime ? timer->fetchNextEvent() : nullptr;
612     constexpr char kernelName[] = "k_pruneonly";
613     const auto     kernel       = plist->haveFreshList ? nbnxn_kernel_prune_cuda<true> : nbnxn_kernel_prune_cuda<false>;
614     const auto     kernelArgs   = prepareGpuKernelArguments(kernel, config, adat, nbp, plist, &numParts, &part);
615     launchGpuKernel(kernel, config, timingEvent, kernelName, kernelArgs);
616
617     /* TODO: consider a more elegant way to track which kernel has been called
618        (combined or separate 1st pass prune, rolling prune). */
619     if (plist->haveFreshList)
620     {
621         plist->haveFreshList                   = false;
622         /* Mark that pruning has been done */
623         nb->timers->interaction[iloc].didPrune = true;
624     }
625     else
626     {
627         /* Mark that rolling pruning has been done */
628         nb->timers->interaction[iloc].didRollingPrune = true;
629     }
630
631     if (bDoTime)
632     {
633         timer->closeTimingRegion(stream);
634     }
635
636     if (GMX_NATIVE_WINDOWS)
637     {
638         /* Windows: force flushing WDDM queue */
639         cudaStreamQuery(stream);
640     }
641 }
642
643 void gpu_launch_cpyback(gmx_nbnxn_cuda_t       *nb,
644                         nbnxn_atomdata_t       *nbatom,
645                         const gmx::ForceFlags  &forceFlags,
646                         const AtomLocality      atomLocality,
647                         const bool              copyBackNbForce)
648 {
649     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
650
651     cudaError_t stat;
652     int         adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
653
654     /* determine interaction locality from atom locality */
655     const InteractionLocality iloc = gpuAtomToInteractionLocality(atomLocality);
656
657     /* extract the data */
658     cu_atomdata_t   *adat    = nb->atdat;
659     cu_timers_t     *t       = nb->timers;
660     bool             bDoTime = nb->bDoTime;
661     cudaStream_t     stream  = nb->stream[iloc];
662
663     /* don't launch non-local copy-back if there was no non-local work to do */
664     if ((iloc == InteractionLocality::NonLocal) && !haveGpuShortRangeWork(*nb, iloc))
665     {
666         return;
667     }
668
669     getGpuAtomRange(adat, atomLocality, &adat_begin, &adat_len);
670
671     /* beginning of timed D2H section */
672     if (bDoTime)
673     {
674         t->xf[atomLocality].nb_d2h.openTimingRegion(stream);
675     }
676
677     /* With DD the local D2H transfer can only start after the non-local
678        kernel has finished. */
679     if (iloc == InteractionLocality::Local && nb->bUseTwoStreams)
680     {
681         stat = cudaStreamWaitEvent(stream, nb->nonlocal_done, 0);
682         CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
683     }
684
685     /* DtoH f */
686     if (copyBackNbForce)
687     {
688         cu_copy_D2H_async(nbatom->out[0].f.data() + adat_begin * 3, adat->f + adat_begin,
689                           (adat_len)*sizeof(*adat->f), stream);
690     }
691
692     /* After the non-local D2H is launched the nonlocal_done event can be
693        recorded which signals that the local D2H can proceed. This event is not
694        placed after the non-local kernel because we want the non-local data
695        back first. */
696     if (iloc == InteractionLocality::NonLocal)
697     {
698         stat = cudaEventRecord(nb->nonlocal_done, stream);
699         CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
700     }
701
702     /* only transfer energies in the local stream */
703     if (iloc == InteractionLocality::Local)
704     {
705         /* DtoH fshift when virial is needed */
706         if (forceFlags.computeVirial)
707         {
708             cu_copy_D2H_async(nb->nbst.fshift, adat->fshift,
709                               SHIFTS * sizeof(*nb->nbst.fshift), stream);
710         }
711
712         /* DtoH energies */
713         if (forceFlags.computeEnergy)
714         {
715             cu_copy_D2H_async(nb->nbst.e_lj, adat->e_lj,
716                               sizeof(*nb->nbst.e_lj), stream);
717             cu_copy_D2H_async(nb->nbst.e_el, adat->e_el,
718                               sizeof(*nb->nbst.e_el), stream);
719         }
720     }
721
722     if (bDoTime)
723     {
724         t->xf[atomLocality].nb_d2h.closeTimingRegion(stream);
725     }
726 }
727
728 void cuda_set_cacheconfig()
729 {
730     cudaError_t stat;
731
732     for (int i = 0; i < eelCuNR; i++)
733     {
734         for (int j = 0; j < evdwCuNR; j++)
735         {
736             /* Default kernel 32/32 kB Shared/L1 */
737             cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferEqual);
738             cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
739             cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferEqual);
740             stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
741             CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");
742         }
743     }
744 }
745
746 /* X buffer operations on GPU: copies coordinates to the GPU in rvec format. */
747 void nbnxn_gpu_copy_x_to_gpu(const Nbnxm::Grid               &grid,
748                              gmx_nbnxn_gpu_t                 *nb,
749                              const Nbnxm::AtomLocality        locality,
750                              const rvec                      *coordinatesHost)
751 {
752     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
753
754     bool                       bDoTime = nb->bDoTime;
755
756     Nbnxm::InteractionLocality interactionLoc            = gpuAtomToInteractionLocality(locality);
757     int                        numCopyAtoms              = grid.srcAtomEnd() - grid.srcAtomBegin();
758     int                        copyAtomStart             = grid.srcAtomBegin();
759
760     cudaStream_t               stream  = nb->stream[interactionLoc];
761
762     // empty domain avoid launching zero-byte copy
763     if (numCopyAtoms == 0)
764     {
765         return;
766     }
767     GMX_ASSERT(coordinatesHost,  "Need a valid host pointer");
768
769     if (bDoTime)
770     {
771         nb->timers->xf[locality].nb_h2d.openTimingRegion(stream);
772     }
773
774     rvec       *devicePtrDest = reinterpret_cast<rvec *> (nb->xrvec[copyAtomStart]);
775     const rvec *devicePtrSrc  = reinterpret_cast<const rvec *> (coordinatesHost[copyAtomStart]);
776     copyToDeviceBuffer(&devicePtrDest, devicePtrSrc, 0, numCopyAtoms,
777                        stream, GpuApiCallBehavior::Async, nullptr);
778
779     if (interactionLoc == Nbnxm::InteractionLocality::Local)
780     {
781         nb->xAvailableOnDevice->markEvent(stream);
782     }
783
784     if (bDoTime)
785     {
786         nb->timers->xf[locality].nb_h2d.closeTimingRegion(stream);
787     }
788 }
789
790 DeviceBuffer<float> nbnxn_gpu_get_x_gpu(gmx_nbnxn_gpu_t *nb)
791 {
792     return reinterpret_cast< DeviceBuffer<float> >(nb->xrvec);
793 }
794
795 /* X buffer operations on GPU: performs conversion from rvec to nb format. */
796 void nbnxn_gpu_x_to_nbat_x(const Nbnxm::Grid               &grid,
797                            bool                             setFillerCoords,
798                            gmx_nbnxn_gpu_t                 *nb,
799                            DeviceBuffer<float>              coordinatesDevice,
800                            const Nbnxm::AtomLocality        locality,
801                            int                              gridId,
802                            int                              numColumnsMax)
803 {
804     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
805
806     cu_atomdata_t             *adat    = nb->atdat;
807
808     const int                  numColumns                = grid.numColumns();
809     const int                  cellOffset                = grid.cellOffset();
810     const int                  numAtomsPerCell           = grid.numAtomsPerCell();
811     Nbnxm::InteractionLocality interactionLoc            = gpuAtomToInteractionLocality(locality);
812
813     cudaStream_t               stream  = nb->stream[interactionLoc];
814
815     int numAtoms = grid.srcAtomEnd() - grid.srcAtomBegin();
816     // avoid empty kernel launch, skip to inserting stream dependency
817     if (numAtoms != 0)
818     {
819         // TODO: This will only work with CUDA
820         GMX_ASSERT(coordinatesDevice, "Need a valid device pointer");
821
822         KernelLaunchConfig config;
823         config.blockSize[0]     = c_bufOpsThreadsPerBlock;
824         config.blockSize[1]     = 1;
825         config.blockSize[2]     = 1;
826         config.gridSize[0]      = (grid.numCellsColumnMax()*numAtomsPerCell + c_bufOpsThreadsPerBlock - 1)/c_bufOpsThreadsPerBlock;
827         config.gridSize[1]      = numColumns;
828         config.gridSize[2]      = 1;
829         GMX_ASSERT(config.gridSize[0] > 0, "Can not have empty grid, early return above avoids this");
830         config.sharedMemorySize = 0;
831         config.stream           = stream;
832
833         auto       kernelFn       = nbnxn_gpu_x_to_nbat_x_kernel;
834         float     *xqPtr          = &(adat->xq->x);
835         const int *d_atomIndices  = nb->atomIndices;
836         const int *d_cxy_na       = &nb->cxy_na[numColumnsMax*gridId];
837         const int *d_cxy_ind      = &nb->cxy_ind[numColumnsMax*gridId];
838         const auto kernelArgs     = prepareGpuKernelArguments(kernelFn, config,
839                                                               &numColumns,
840                                                               &xqPtr,
841                                                               &setFillerCoords,
842                                                               &coordinatesDevice,
843                                                               &d_atomIndices,
844                                                               &d_cxy_na,
845                                                               &d_cxy_ind,
846                                                               &cellOffset,
847                                                               &numAtomsPerCell);
848         launchGpuKernel(kernelFn, config, nullptr, "XbufferOps", kernelArgs);
849     }
850
851     // TODO: note that this is not necessary when there are no local atoms, that is:
852     // (numAtoms == 0 && interactionLoc == InteractionLocality::Local)
853     // but for now we avoid that optimization
854     nbnxnInsertNonlocalGpuDependency(nb, interactionLoc);
855 }
856
857 /* F buffer operations on GPU: performs force summations and conversion from nb to rvec format. */
858 void nbnxn_gpu_add_nbat_f_to_f(const AtomLocality               atomLocality,
859                                DeviceBuffer<float>              totalForcesDevice,
860                                gmx_nbnxn_gpu_t                 *nb,
861                                void                            *pmeForcesDevice,
862                                GpuEventSynchronizer            *pmeForcesReady,
863                                int                              atomStart,
864                                int                              numAtoms,
865                                bool                             useGpuFPmeReduction,
866                                bool                             accumulateForce)
867 {
868     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
869     GMX_ASSERT(numAtoms != 0, "Cannot call function with no atoms");
870     GMX_ASSERT(totalForcesDevice, "Need a valid totalForcesDevice pointer");
871
872     const InteractionLocality iLocality     = gpuAtomToInteractionLocality(atomLocality);
873     cudaStream_t              stream        = nb->stream[iLocality];
874     cu_atomdata_t            *adat          = nb->atdat;
875
876     if (useGpuFPmeReduction)
877     {
878         //Stream must wait for PME force completion
879         pmeForcesReady->enqueueWaitEvent(stream);
880     }
881
882     /* launch kernel */
883
884     KernelLaunchConfig config;
885     config.blockSize[0]     = c_bufOpsThreadsPerBlock;
886     config.blockSize[1]     = 1;
887     config.blockSize[2]     = 1;
888     config.gridSize[0]      = ((numAtoms+1)+c_bufOpsThreadsPerBlock-1)/c_bufOpsThreadsPerBlock;
889     config.gridSize[1]      = 1;
890     config.gridSize[2]      = 1;
891     config.sharedMemorySize = 0;
892     config.stream           = stream;
893
894     auto  kernelFn = accumulateForce ?
895         nbnxn_gpu_add_nbat_f_to_f_kernel<true, false> :
896         nbnxn_gpu_add_nbat_f_to_f_kernel<false, false>;
897
898     if (useGpuFPmeReduction)
899     {
900         GMX_ASSERT(pmeForcesDevice, "Need a valid pmeForcesDevice pointer");
901         kernelFn = accumulateForce ?
902             nbnxn_gpu_add_nbat_f_to_f_kernel<true, true> :
903             nbnxn_gpu_add_nbat_f_to_f_kernel<false, true>;
904     }
905
906     const float3     *d_fNB    = adat->f;
907     const float3     *d_fPme   = (float3*) pmeForcesDevice;
908     float3           *d_fTotal = (float3*) totalForcesDevice;
909     const int        *d_cell   = nb->cell;
910
911     const auto        kernelArgs   = prepareGpuKernelArguments(kernelFn, config,
912                                                                &d_fNB,
913                                                                &d_fPme,
914                                                                &d_fTotal,
915                                                                &d_cell,
916                                                                &atomStart,
917                                                                &numAtoms);
918
919     launchGpuKernel(kernelFn, config, nullptr, "FbufferOps", kernelArgs);
920
921 }
922
923 DeviceBuffer<float> nbnxn_gpu_get_f_gpu(gmx_nbnxn_gpu_t *nb)
924 {
925     return reinterpret_cast< DeviceBuffer<float> >(nb->frvec);
926 }
927
928 void nbnxn_launch_copy_f_to_gpu(const AtomLocality               atomLocality,
929                                 const Nbnxm::GridSet            &gridSet,
930                                 gmx_nbnxn_gpu_t                 *nb,
931                                 rvec                            *f)
932 {
933     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
934
935     const InteractionLocality iLocality = gpuAtomToInteractionLocality(atomLocality);
936     cudaStream_t              stream    = nb->stream[iLocality];
937
938     bool                      bDoTime = nb->bDoTime;
939     cu_timers_t              *t       = nb->timers;
940
941     int                       atomStart = 0, numCopyAtoms = 0;
942
943     nbnxn_get_atom_range(atomLocality, gridSet, &atomStart, &numCopyAtoms);
944
945     // Avoiding launching copy with no work
946     if (numCopyAtoms == 0)
947     {
948         return;
949     }
950     GMX_ASSERT(f, "Need a valid f pointer");
951
952     if (bDoTime)
953     {
954         t->xf[atomLocality].nb_h2d.openTimingRegion(stream);
955     }
956
957     rvec       *ptrDest  = reinterpret_cast<rvec *> (nb->frvec[atomStart]);
958     rvec       *ptrSrc   = reinterpret_cast<rvec *> (f[atomStart]);
959     //copyToDeviceBuffer(&ptrDest, ptrSrc, 0, numCopyAtoms,
960     //                   stream, GpuApiCallBehavior::Async, nullptr);
961     //TODO use above API call rather than direct memcpy when force has been implemented in a hostvector
962     cudaMemcpyAsync(ptrDest, ptrSrc, numCopyAtoms*sizeof(rvec), cudaMemcpyHostToDevice,
963                     stream);
964
965     if (bDoTime)
966     {
967         t->xf[atomLocality].nb_h2d.closeTimingRegion(stream);
968     }
969
970     return;
971 }
972
973 void nbnxn_launch_copy_f_from_gpu(const AtomLocality               atomLocality,
974                                   const Nbnxm::GridSet            &gridSet,
975                                   gmx_nbnxn_gpu_t                 *nb,
976                                   rvec                            *f)
977 {
978     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
979
980     const InteractionLocality iLocality = gpuAtomToInteractionLocality(atomLocality);
981     cudaStream_t              stream    = nb->stream[iLocality];
982
983     bool                      bDoTime = nb->bDoTime;
984     cu_timers_t              *t       = nb->timers;
985     int                       atomStart, numCopyAtoms;
986
987     nbnxn_get_atom_range(atomLocality, gridSet, &atomStart, &numCopyAtoms);
988
989     // Avoiding launching copy with no work
990     if (numCopyAtoms == 0)
991     {
992         return;
993     }
994     GMX_ASSERT(f, "Need a valid f pointer");
995
996     if (bDoTime)
997     {
998         t->xf[atomLocality].nb_d2h.openTimingRegion(stream);
999     }
1000
1001     GMX_ASSERT(nb->frvec,  "Need a valid nb->frvec pointer");
1002     rvec       *ptrDest = reinterpret_cast<rvec *> (f[atomStart]);
1003     rvec       *ptrSrc  = reinterpret_cast<rvec *> (nb->frvec[atomStart]);
1004     //copyFromDeviceBuffer(ptrDest, &ptrSrc, 0, numCopyAtoms,
1005     //                   stream, GpuApiCallBehavior::Async, nullptr);
1006     //TODO use above API call rather than direct memcpy when force has been implemented in a hostvector
1007     cudaMemcpyAsync(ptrDest, ptrSrc, numCopyAtoms*sizeof(rvec), cudaMemcpyDeviceToHost,
1008                     stream);
1009
1010     if (bDoTime)
1011     {
1012         t->xf[atomLocality].nb_d2h.closeTimingRegion(stream);
1013     }
1014
1015     return;
1016 }
1017
1018 void nbnxn_launch_copy_x_from_gpu(const AtomLocality               atomLocality,
1019                                   const Nbnxm::GridSet            &gridSet,
1020                                   gmx_nbnxn_gpu_t                 *nb,
1021                                   rvec                            *x)
1022 {
1023     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
1024     GMX_ASSERT(x,  "Need a valid x pointer");
1025
1026     const InteractionLocality iLocality = gpuAtomToInteractionLocality(atomLocality);
1027     cudaStream_t              stream    = nb->stream[iLocality];
1028
1029     bool                      bDoTime = nb->bDoTime;
1030     cu_timers_t              *t       = nb->timers;
1031     int                       atomStart, nAtoms;
1032
1033     nbnxn_get_atom_range(atomLocality, gridSet, &atomStart, &nAtoms);
1034
1035     if (bDoTime)
1036     {
1037         t->xf[atomLocality].nb_d2h.openTimingRegion(stream);
1038     }
1039
1040     GMX_ASSERT(nb->xrvec,  "Need a valid nb->xrvec pointer");
1041     rvec       *ptrDest = reinterpret_cast<rvec *> (x[atomStart]);
1042     rvec       *ptrSrc  = reinterpret_cast<rvec *> (nb->xrvec[atomStart]);
1043     copyFromDeviceBuffer(ptrDest, &ptrSrc, 0, nAtoms,
1044                          stream, GpuApiCallBehavior::Async, stream);
1045
1046     if (atomLocality == AtomLocality::NonLocal)
1047     {
1048         nb->xNonLocalCopyD2HDone->markEvent(stream);
1049     }
1050
1051     if (bDoTime)
1052     {
1053         t->xf[atomLocality].nb_d2h.closeTimingRegion(stream);
1054     }
1055
1056     return;
1057 }
1058
1059 void nbnxn_wait_for_gpu_force_reduction(const AtomLocality      gmx_unused atomLocality,
1060                                         gmx_nbnxn_gpu_t                   *nb)
1061 {
1062     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
1063
1064     const InteractionLocality iLocality = gpuAtomToInteractionLocality(atomLocality);
1065
1066     cudaStream_t              stream    = nb->stream[iLocality];
1067
1068     cudaStreamSynchronize(stream);
1069
1070 }
1071
1072 void* nbnxn_get_gpu_xrvec(gmx_nbnxn_gpu_t *gpu_nbv)
1073 {
1074     return static_cast<void *> (gpu_nbv->xrvec);
1075 }
1076
1077 void* nbnxn_get_gpu_frvec(gmx_nbnxn_gpu_t *gpu_nbv)
1078 {
1079     return static_cast<void *> (gpu_nbv->frvec);
1080 }
1081
1082 void* nbnxn_get_x_on_device_event(const gmx_nbnxn_cuda_t   *nb)
1083 {
1084     return static_cast<void*> (nb->xAvailableOnDevice);
1085 }
1086
1087 void nbnxn_wait_nonlocal_x_copy_D2H_done(gmx_nbnxn_cuda_t   *nb)
1088 {
1089     nb->xNonLocalCopyD2HDone->waitForEvent();
1090 }
1091
1092 void nbnxn_stream_local_wait_for_nonlocal(gmx_nbnxn_cuda_t   *nb)
1093 {
1094     cudaStream_t         localStream     = nb->stream[InteractionLocality::Local];
1095     cudaStream_t         nonLocalStream  = nb->stream[InteractionLocality::NonLocal];
1096
1097     GpuEventSynchronizer event;
1098     event.markEvent(nonLocalStream);
1099     event.enqueueWaitEvent(localStream);
1100 }
1101
1102 } // namespace Nbnxm