Merge branch release-2019 into release-2020
[alexxy/gromacs.git] / src / gromacs / nbnxm / opencl / nbnxm_ocl.cpp
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 /*! \internal \file
36  *  \brief Define OpenCL implementation of nbnxm_gpu.h
37  *
38  *  \author Anca Hamuraru <anca@streamcomputing.eu>
39  *  \author Teemu Virolainen <teemu@streamcomputing.eu>
40  *  \author Dimitrios Karkoulis <dimitris.karkoulis@gmail.com>
41  *  \author Szilárd Páll <pall.szilard@gmail.com>
42  *  \ingroup module_nbnxm
43  *
44  *  TODO (psz):
45  *  - Add a static const cl_uint c_pruneKernelWorkDim / c_nbnxnKernelWorkDim = 3;
46  *  - Rework the copying of OCL data structures done before every invocation of both
47  *    nb and prune kernels (using fillin_ocl_structures); also consider at the same
48  *    time calling clSetKernelArg only on the updated parameters (if tracking changed
49  *    parameters is feasible);
50  *  - Consider using the event_wait_list argument to clEnqueueNDRangeKernel to mark
51  *    dependencies on the kernel launched: e.g. the non-local nb kernel's dependency
52  *    on the misc_ops_and_local_H2D_done event could be better expressed this way.
53  *
54  *  - Consider extracting common sections of the OpenCL and CUDA nbnxn logic, e.g:
55  *    - in nbnxn_gpu_launch_kernel_pruneonly() the pre- and post-kernel launch logic
56  *      is identical in the two implementations, so a 3-way split might allow sharing
57  *      code;
58  *    -
59  *
60  */
61 #include "gmxpre.h"
62
63 #include <assert.h>
64 #include <stdlib.h>
65
66 #if defined(_MSVC)
67 #    include <limits>
68 #endif
69
70 #include "thread_mpi/atomic.h"
71
72 #include "gromacs/gpu_utils/gputraits_ocl.h"
73 #include "gromacs/gpu_utils/oclutils.h"
74 #include "gromacs/hardware/hw_info.h"
75 #include "gromacs/mdtypes/simulation_workload.h"
76 #include "gromacs/nbnxm/atomdata.h"
77 #include "gromacs/nbnxm/gpu_common.h"
78 #include "gromacs/nbnxm/gpu_common_utils.h"
79 #include "gromacs/nbnxm/gpu_data_mgmt.h"
80 #include "gromacs/nbnxm/nbnxm.h"
81 #include "gromacs/nbnxm/nbnxm_gpu.h"
82 #include "gromacs/nbnxm/pairlist.h"
83 #include "gromacs/pbcutil/ishift.h"
84 #include "gromacs/timing/gpu_timing.h"
85 #include "gromacs/utility/cstringutil.h"
86 #include "gromacs/utility/fatalerror.h"
87 #include "gromacs/utility/gmxassert.h"
88
89 #include "nbnxm_ocl_internal.h"
90 #include "nbnxm_ocl_types.h"
91
92 namespace Nbnxm
93 {
94
95 /*! \brief Convenience constants */
96 //@{
97 static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
98 static const int c_clSize          = c_nbnxnGpuClusterSize;
99 //@}
100
101
102 /*! \brief Validates the input global work size parameter.
103  */
104 static inline void validate_global_work_size(const KernelLaunchConfig& config,
105                                              int                       work_dim,
106                                              const gmx_device_info_t*  dinfo)
107 {
108     cl_uint device_size_t_size_bits;
109     cl_uint host_size_t_size_bits;
110
111     assert(dinfo);
112
113     size_t global_work_size[3];
114     GMX_ASSERT(work_dim <= 3, "Not supporting hyper-grids just yet");
115     for (int i = 0; i < work_dim; i++)
116     {
117         global_work_size[i] = config.blockSize[i] * config.gridSize[i];
118     }
119
120     /* Each component of a global_work_size must not exceed the range given by the
121        sizeof(device size_t) for the device on which the kernel execution will
122        be enqueued. See:
123        https://www.khronos.org/registry/cl/sdk/1.0/docs/man/xhtml/clEnqueueNDRangeKernel.html
124      */
125     device_size_t_size_bits = dinfo->adress_bits;
126     host_size_t_size_bits   = static_cast<cl_uint>(sizeof(size_t) * 8);
127
128     /* If sizeof(host size_t) <= sizeof(device size_t)
129             => global_work_size components will always be valid
130        else
131             => get device limit for global work size and
132             compare it against each component of global_work_size.
133      */
134     if (host_size_t_size_bits > device_size_t_size_bits)
135     {
136         size_t device_limit;
137
138         device_limit = (1ULL << device_size_t_size_bits) - 1;
139
140         for (int i = 0; i < work_dim; i++)
141         {
142             if (global_work_size[i] > device_limit)
143             {
144                 gmx_fatal(
145                         FARGS,
146                         "Watch out, the input system is too large to simulate!\n"
147                         "The number of nonbonded work units (=number of super-clusters) exceeds the"
148                         "device capabilities. Global work size limit exceeded (%zu > %zu)!",
149                         global_work_size[i], device_limit);
150             }
151         }
152     }
153 }
154
155 /* Constant arrays listing non-bonded kernel function names. The arrays are
156  * organized in 2-dim arrays by: electrostatics and VDW type.
157  *
158  *  Note that the row- and column-order of function pointers has to match the
159  *  order of corresponding enumerated electrostatics and vdw types, resp.,
160  *  defined in nbnxm_ocl_types.h.
161  */
162
163 /*! \brief Force-only kernel function names. */
164 static const char* nb_kfunc_noener_noprune_ptr[eelOclNR][evdwOclNR] = {
165     { "nbnxn_kernel_ElecCut_VdwLJ_F_opencl", "nbnxn_kernel_ElecCut_VdwLJCombGeom_F_opencl",
166       "nbnxn_kernel_ElecCut_VdwLJCombLB_F_opencl", "nbnxn_kernel_ElecCut_VdwLJFsw_F_opencl",
167       "nbnxn_kernel_ElecCut_VdwLJPsw_F_opencl", "nbnxn_kernel_ElecCut_VdwLJEwCombGeom_F_opencl",
168       "nbnxn_kernel_ElecCut_VdwLJEwCombLB_F_opencl" },
169     { "nbnxn_kernel_ElecRF_VdwLJ_F_opencl", "nbnxn_kernel_ElecRF_VdwLJCombGeom_F_opencl",
170       "nbnxn_kernel_ElecRF_VdwLJCombLB_F_opencl", "nbnxn_kernel_ElecRF_VdwLJFsw_F_opencl",
171       "nbnxn_kernel_ElecRF_VdwLJPsw_F_opencl", "nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_opencl",
172       "nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_opencl" },
173     { "nbnxn_kernel_ElecEwQSTab_VdwLJ_F_opencl", "nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_F_opencl",
174       "nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_F_opencl", "nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_opencl",
175       "nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_opencl",
176       "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_F_opencl",
177       "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_F_opencl" },
178     { "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_opencl",
179       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_F_opencl",
180       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_F_opencl",
181       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_opencl",
182       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_opencl",
183       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_F_opencl",
184       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_F_opencl" },
185     { "nbnxn_kernel_ElecEw_VdwLJ_F_opencl", "nbnxn_kernel_ElecEw_VdwLJCombGeom_F_opencl",
186       "nbnxn_kernel_ElecEw_VdwLJCombLB_F_opencl", "nbnxn_kernel_ElecEw_VdwLJFsw_F_opencl",
187       "nbnxn_kernel_ElecEw_VdwLJPsw_F_opencl", "nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_opencl",
188       "nbnxn_kernel_ElecEw_VdwLJEwCombLB_F_opencl" },
189     { "nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_opencl",
190       "nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_opencl",
191       "nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_opencl",
192       "nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_opencl", "nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_opencl",
193       "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_opencl",
194       "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_F_opencl" }
195 };
196
197 /*! \brief Force + energy kernel function pointers. */
198 static const char* nb_kfunc_ener_noprune_ptr[eelOclNR][evdwOclNR] = {
199     { "nbnxn_kernel_ElecCut_VdwLJ_VF_opencl", "nbnxn_kernel_ElecCut_VdwLJCombGeom_VF_opencl",
200       "nbnxn_kernel_ElecCut_VdwLJCombLB_VF_opencl", "nbnxn_kernel_ElecCut_VdwLJFsw_VF_opencl",
201       "nbnxn_kernel_ElecCut_VdwLJPsw_VF_opencl", "nbnxn_kernel_ElecCut_VdwLJEwCombGeom_VF_opencl",
202       "nbnxn_kernel_ElecCut_VdwLJEwCombLB_VF_opencl" },
203     { "nbnxn_kernel_ElecRF_VdwLJ_VF_opencl", "nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_opencl",
204       "nbnxn_kernel_ElecRF_VdwLJCombLB_VF_opencl", "nbnxn_kernel_ElecRF_VdwLJFsw_VF_opencl",
205       "nbnxn_kernel_ElecRF_VdwLJPsw_VF_opencl", "nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_opencl",
206       "nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_opencl" },
207     { "nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_opencl", "nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_VF_opencl",
208       "nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_VF_opencl",
209       "nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_opencl", "nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_opencl",
210       "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_VF_opencl",
211       "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_VF_opencl" },
212     { "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_opencl",
213       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_VF_opencl",
214       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_VF_opencl",
215       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_opencl",
216       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_opencl",
217       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_VF_opencl",
218       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_VF_opencl" },
219     { "nbnxn_kernel_ElecEw_VdwLJ_VF_opencl", "nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_opencl",
220       "nbnxn_kernel_ElecEw_VdwLJCombLB_VF_opencl", "nbnxn_kernel_ElecEw_VdwLJFsw_VF_opencl",
221       "nbnxn_kernel_ElecEw_VdwLJPsw_VF_opencl", "nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_opencl",
222       "nbnxn_kernel_ElecEw_VdwLJEwCombLB_VF_opencl" },
223     { "nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_opencl",
224       "nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_opencl",
225       "nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_opencl",
226       "nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_opencl",
227       "nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_opencl",
228       "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_opencl",
229       "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_VF_opencl" }
230 };
231
232 /*! \brief Force + pruning kernel function pointers. */
233 static const char* nb_kfunc_noener_prune_ptr[eelOclNR][evdwOclNR] = {
234     { "nbnxn_kernel_ElecCut_VdwLJ_F_prune_opencl",
235       "nbnxn_kernel_ElecCut_VdwLJCombGeom_F_prune_opencl",
236       "nbnxn_kernel_ElecCut_VdwLJCombLB_F_prune_opencl",
237       "nbnxn_kernel_ElecCut_VdwLJFsw_F_prune_opencl", "nbnxn_kernel_ElecCut_VdwLJPsw_F_prune_opencl",
238       "nbnxn_kernel_ElecCut_VdwLJEwCombGeom_F_prune_opencl",
239       "nbnxn_kernel_ElecCut_VdwLJEwCombLB_F_prune_opencl" },
240     { "nbnxn_kernel_ElecRF_VdwLJ_F_prune_opencl", "nbnxn_kernel_ElecRF_VdwLJCombGeom_F_prune_opencl",
241       "nbnxn_kernel_ElecRF_VdwLJCombLB_F_prune_opencl",
242       "nbnxn_kernel_ElecRF_VdwLJFsw_F_prune_opencl", "nbnxn_kernel_ElecRF_VdwLJPsw_F_prune_opencl",
243       "nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_prune_opencl",
244       "nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_prune_opencl" },
245     { "nbnxn_kernel_ElecEwQSTab_VdwLJ_F_prune_opencl",
246       "nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_F_prune_opencl",
247       "nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_F_prune_opencl",
248       "nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_prune_opencl",
249       "nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_prune_opencl",
250       "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_F_prune_opencl",
251       "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_F_prune_opencl" },
252     { "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_prune_opencl",
253       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_F_prune_opencl",
254       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_F_prune_opencl",
255       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_prune_opencl",
256       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_prune_opencl",
257       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_F_prune_opencl",
258       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_F_prune_opencl" },
259     { "nbnxn_kernel_ElecEw_VdwLJ_F_prune_opencl", "nbnxn_kernel_ElecEw_VdwLJCombGeom_F_prune_opencl",
260       "nbnxn_kernel_ElecEw_VdwLJCombLB_F_prune_opencl",
261       "nbnxn_kernel_ElecEw_VdwLJFsw_F_prune_opencl", "nbnxn_kernel_ElecEw_VdwLJPsw_F_prune_opencl",
262       "nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_prune_opencl",
263       "nbnxn_kernel_ElecEw_VdwLJEwCombLB_F_prune_opencl" },
264     { "nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_prune_opencl",
265       "nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_prune_opencl",
266       "nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_prune_opencl",
267       "nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_prune_opencl",
268       "nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_prune_opencl",
269       "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_prune_opencl",
270       "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_F_prune_opencl" }
271 };
272
273 /*! \brief Force + energy + pruning kernel function pointers. */
274 static const char* nb_kfunc_ener_prune_ptr[eelOclNR][evdwOclNR] = {
275     { "nbnxn_kernel_ElecCut_VdwLJ_VF_prune_opencl",
276       "nbnxn_kernel_ElecCut_VdwLJCombGeom_VF_prune_opencl",
277       "nbnxn_kernel_ElecCut_VdwLJCombLB_VF_prune_opencl",
278       "nbnxn_kernel_ElecCut_VdwLJFsw_VF_prune_opencl",
279       "nbnxn_kernel_ElecCut_VdwLJPsw_VF_prune_opencl",
280       "nbnxn_kernel_ElecCut_VdwLJEwCombGeom_VF_prune_opencl",
281       "nbnxn_kernel_ElecCut_VdwLJEwCombLB_VF_prune_opencl" },
282     { "nbnxn_kernel_ElecRF_VdwLJ_VF_prune_opencl",
283       "nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_prune_opencl",
284       "nbnxn_kernel_ElecRF_VdwLJCombLB_VF_prune_opencl",
285       "nbnxn_kernel_ElecRF_VdwLJFsw_VF_prune_opencl", "nbnxn_kernel_ElecRF_VdwLJPsw_VF_prune_opencl",
286       "nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_prune_opencl",
287       "nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_prune_opencl" },
288     { "nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_prune_opencl",
289       "nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_VF_prune_opencl",
290       "nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_VF_prune_opencl",
291       "nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_prune_opencl",
292       "nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_prune_opencl",
293       "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_VF_prune_opencl",
294       "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_VF_prune_opencl" },
295     { "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_prune_opencl",
296       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_VF_prune_opencl",
297       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_VF_prune_opencl",
298       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_prune_opencl",
299       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_prune_opencl",
300       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_VF_prune_opencl",
301       "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_VF_prune_opencl" },
302     { "nbnxn_kernel_ElecEw_VdwLJ_VF_prune_opencl",
303       "nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_prune_opencl",
304       "nbnxn_kernel_ElecEw_VdwLJCombLB_VF_prune_opencl",
305       "nbnxn_kernel_ElecEw_VdwLJFsw_VF_prune_opencl", "nbnxn_kernel_ElecEw_VdwLJPsw_VF_prune_opencl",
306       "nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_prune_opencl",
307       "nbnxn_kernel_ElecEw_VdwLJEwCombLB_VF_prune_opencl" },
308     { "nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_prune_opencl",
309       "nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_prune_opencl",
310       "nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_prune_opencl",
311       "nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_prune_opencl",
312       "nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_prune_opencl",
313       "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_prune_opencl",
314       "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_VF_prune_opencl" }
315 };
316
317 /*! \brief Return a pointer to the prune kernel version to be executed at the current invocation.
318  *
319  * \param[in] kernel_pruneonly  array of prune kernel objects
320  * \param[in] firstPrunePass    true if the first pruning pass is being executed
321  */
322 static inline cl_kernel selectPruneKernel(cl_kernel kernel_pruneonly[], bool firstPrunePass)
323 {
324     cl_kernel* kernelPtr;
325
326     if (firstPrunePass)
327     {
328         kernelPtr = &(kernel_pruneonly[epruneFirst]);
329     }
330     else
331     {
332         kernelPtr = &(kernel_pruneonly[epruneRolling]);
333     }
334     // TODO: consider creating the prune kernel object here to avoid a
335     // clCreateKernel for the rolling prune kernel if this is not needed.
336     return *kernelPtr;
337 }
338
339 /*! \brief Return a pointer to the kernel version to be executed at the current step.
340  *  OpenCL kernel objects are cached in nb. If the requested kernel is not
341  *  found in the cache, it will be created and the cache will be updated.
342  */
343 static inline cl_kernel select_nbnxn_kernel(gmx_nbnxn_ocl_t* nb, int eeltype, int evdwtype, bool bDoEne, bool bDoPrune)
344 {
345     const char* kernel_name_to_run;
346     cl_kernel*  kernel_ptr;
347     cl_int      cl_error;
348
349     assert(eeltype < eelOclNR);
350     assert(evdwtype < evdwOclNR);
351
352     if (bDoEne)
353     {
354         if (bDoPrune)
355         {
356             kernel_name_to_run = nb_kfunc_ener_prune_ptr[eeltype][evdwtype];
357             kernel_ptr         = &(nb->kernel_ener_prune_ptr[eeltype][evdwtype]);
358         }
359         else
360         {
361             kernel_name_to_run = nb_kfunc_ener_noprune_ptr[eeltype][evdwtype];
362             kernel_ptr         = &(nb->kernel_ener_noprune_ptr[eeltype][evdwtype]);
363         }
364     }
365     else
366     {
367         if (bDoPrune)
368         {
369             kernel_name_to_run = nb_kfunc_noener_prune_ptr[eeltype][evdwtype];
370             kernel_ptr         = &(nb->kernel_noener_prune_ptr[eeltype][evdwtype]);
371         }
372         else
373         {
374             kernel_name_to_run = nb_kfunc_noener_noprune_ptr[eeltype][evdwtype];
375             kernel_ptr         = &(nb->kernel_noener_noprune_ptr[eeltype][evdwtype]);
376         }
377     }
378
379     if (nullptr == kernel_ptr[0])
380     {
381         *kernel_ptr = clCreateKernel(nb->dev_rundata->program, kernel_name_to_run, &cl_error);
382         assert(cl_error == CL_SUCCESS);
383     }
384     // TODO: handle errors
385
386     return *kernel_ptr;
387 }
388
389 /*! \brief Calculates the amount of shared memory required by the nonbonded kernel in use.
390  */
391 static inline int calc_shmem_required_nonbonded(int vdwType, bool bPrefetchLjParam)
392 {
393     int shmem;
394
395     /* size of shmem (force-buffers/xq/atom type preloading) */
396     /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
397     /* i-atom x+q in shared memory */
398     shmem = c_numClPerSupercl * c_clSize * sizeof(float) * 4; /* xqib */
399     /* cj in shared memory, for both warps separately
400      * TODO: in the "nowarp kernels we load cj only once  so the factor 2 is not needed.
401      */
402     shmem += 2 * c_nbnxnGpuJgroupSize * sizeof(int); /* cjs  */
403     if (bPrefetchLjParam)
404     {
405         if (useLjCombRule(vdwType))
406         {
407             /* i-atom LJ combination parameters in shared memory */
408             shmem += c_numClPerSupercl * c_clSize * 2 * sizeof(float); /* atib abused for ljcp, float2 */
409         }
410         else
411         {
412             /* i-atom types in shared memory */
413             shmem += c_numClPerSupercl * c_clSize * sizeof(int); /* atib */
414         }
415     }
416     /* force reduction buffers in shared memory */
417     shmem += c_clSize * c_clSize * 3 * sizeof(float); /* f_buf */
418     /* Warp vote. In fact it must be * number of warps in block.. */
419     shmem += sizeof(cl_uint) * 2; /* warp_any */
420     return shmem;
421 }
422
423 /*! \brief Initializes data structures that are going to be sent to the OpenCL device.
424  *
425  *  The device can't use the same data structures as the host for two main reasons:
426  *  - OpenCL restrictions (pointers are not accepted inside data structures)
427  *  - some host side fields are not needed for the OpenCL kernels.
428  *
429  *  This function is called before the launch of both nbnxn and prune kernels.
430  */
431 static void fillin_ocl_structures(cl_nbparam_t* nbp, cl_nbparam_params_t* nbparams_params)
432 {
433     nbparams_params->coulomb_tab_scale = nbp->coulomb_tab_scale;
434     nbparams_params->c_rf              = nbp->c_rf;
435     nbparams_params->dispersion_shift  = nbp->dispersion_shift;
436     nbparams_params->eeltype           = nbp->eeltype;
437     nbparams_params->epsfac            = nbp->epsfac;
438     nbparams_params->ewaldcoeff_lj     = nbp->ewaldcoeff_lj;
439     nbparams_params->ewald_beta        = nbp->ewald_beta;
440     nbparams_params->rcoulomb_sq       = nbp->rcoulomb_sq;
441     nbparams_params->repulsion_shift   = nbp->repulsion_shift;
442     nbparams_params->rlistOuter_sq     = nbp->rlistOuter_sq;
443     nbparams_params->rvdw_sq           = nbp->rvdw_sq;
444     nbparams_params->rlistInner_sq     = nbp->rlistInner_sq;
445     nbparams_params->rvdw_switch       = nbp->rvdw_switch;
446     nbparams_params->sh_ewald          = nbp->sh_ewald;
447     nbparams_params->sh_lj_ewald       = nbp->sh_lj_ewald;
448     nbparams_params->two_k_rf          = nbp->two_k_rf;
449     nbparams_params->vdwtype           = nbp->vdwtype;
450     nbparams_params->vdw_switch        = nbp->vdw_switch;
451 }
452
453 /*! \brief Enqueues a wait for event completion.
454  *
455  * Then it releases the event and sets it to 0.
456  * Don't use this function when more than one wait will be issued for the event.
457  * Equivalent to Cuda Stream Sync. */
458 static void sync_ocl_event(cl_command_queue stream, cl_event* ocl_event)
459 {
460     cl_int gmx_unused cl_error;
461
462     /* Enqueue wait */
463     cl_error = clEnqueueBarrierWithWaitList(stream, 1, ocl_event, nullptr);
464     GMX_RELEASE_ASSERT(CL_SUCCESS == cl_error, ocl_get_error_string(cl_error).c_str());
465
466     /* Release event and reset it to 0. It is ok to release it as enqueuewaitforevents performs implicit retain for events. */
467     cl_error = clReleaseEvent(*ocl_event);
468     assert(CL_SUCCESS == cl_error);
469     *ocl_event = nullptr;
470 }
471
472 /*! \brief Launch asynchronously the xq buffer host to device copy. */
473 void gpu_copy_xq_to_gpu(gmx_nbnxn_ocl_t* nb, const nbnxn_atomdata_t* nbatom, const AtomLocality atomLocality)
474 {
475     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
476
477     const InteractionLocality iloc = gpuAtomToInteractionLocality(atomLocality);
478
479     /* local/nonlocal offset and length used for xq and f */
480     int adat_begin, adat_len;
481
482     cl_atomdata_t*   adat   = nb->atdat;
483     cl_plist_t*      plist  = nb->plist[iloc];
484     cl_timers_t*     t      = nb->timers;
485     cl_command_queue stream = nb->stream[iloc];
486
487     bool bDoTime = (nb->bDoTime) != 0;
488
489     /* Don't launch the non-local H2D copy if there is no dependent
490        work to do: neither non-local nor other (e.g. bonded) work
491        to do that has as input the nbnxn coordinates.
492        Doing the same for the local kernel is more complicated, since the
493        local part of the force array also depends on the non-local kernel.
494        So to avoid complicating the code and to reduce the risk of bugs,
495        we always call the local local x+q copy (and the rest of the local
496        work in nbnxn_gpu_launch_kernel().
497      */
498     if ((iloc == InteractionLocality::NonLocal) && !haveGpuShortRangeWork(*nb, iloc))
499     {
500         plist->haveFreshList = false;
501
502         return;
503     }
504
505     /* calculate the atom data index range based on locality */
506     if (atomLocality == AtomLocality::Local)
507     {
508         adat_begin = 0;
509         adat_len   = adat->natoms_local;
510     }
511     else
512     {
513         adat_begin = adat->natoms_local;
514         adat_len   = adat->natoms - adat->natoms_local;
515     }
516
517     /* beginning of timed HtoD section */
518     if (bDoTime)
519     {
520         t->xf[atomLocality].nb_h2d.openTimingRegion(stream);
521     }
522
523     /* HtoD x, q */
524     ocl_copy_H2D_async(adat->xq, nbatom->x().data() + adat_begin * 4,
525                        adat_begin * sizeof(float) * 4, adat_len * sizeof(float) * 4, stream,
526                        bDoTime ? t->xf[atomLocality].nb_h2d.fetchNextEvent() : nullptr);
527
528     if (bDoTime)
529     {
530         t->xf[atomLocality].nb_h2d.closeTimingRegion(stream);
531     }
532
533     /* When we get here all misc operations issues in the local stream as well as
534        the local xq H2D are done,
535        so we record that in the local stream and wait for it in the nonlocal one. */
536     if (nb->bUseTwoStreams)
537     {
538         if (iloc == InteractionLocality::Local)
539         {
540             cl_int gmx_used_in_debug cl_error = clEnqueueMarkerWithWaitList(
541                     stream, 0, nullptr, &(nb->misc_ops_and_local_H2D_done));
542             assert(CL_SUCCESS == cl_error);
543
544             /* Based on the v1.2 section 5.13 of the OpenCL spec, a flush is needed
545              * in the local stream in order to be able to sync with the above event
546              * from the non-local stream.
547              */
548             cl_error = clFlush(stream);
549             assert(CL_SUCCESS == cl_error);
550         }
551         else
552         {
553             sync_ocl_event(stream, &(nb->misc_ops_and_local_H2D_done));
554         }
555     }
556 }
557
558
559 /*! \brief Launch GPU kernel
560
561    As we execute nonbonded workload in separate queues, before launching
562    the kernel we need to make sure that he following operations have completed:
563    - atomdata allocation and related H2D transfers (every nstlist step);
564    - pair list H2D transfer (every nstlist step);
565    - shift vector H2D transfer (every nstlist step);
566    - force (+shift force and energy) output clearing (every step).
567
568    These operations are issued in the local queue at the beginning of the step
569    and therefore always complete before the local kernel launch. The non-local
570    kernel is launched after the local on the same device/context, so this is
571    inherently scheduled after the operations in the local stream (including the
572    above "misc_ops").
573    However, for the sake of having a future-proof implementation, we use the
574    misc_ops_done event to record the point in time when the above  operations
575    are finished and synchronize with this event in the non-local stream.
576  */
577 void gpu_launch_kernel(gmx_nbnxn_ocl_t* nb, const gmx::StepWorkload& stepWork, const Nbnxm::InteractionLocality iloc)
578 {
579     cl_atomdata_t*   adat   = nb->atdat;
580     cl_nbparam_t*    nbp    = nb->nbparam;
581     cl_plist_t*      plist  = nb->plist[iloc];
582     cl_timers_t*     t      = nb->timers;
583     cl_command_queue stream = nb->stream[iloc];
584
585     bool bDoTime = (nb->bDoTime) != 0;
586
587     cl_nbparam_params_t nbparams_params;
588
589     /* Don't launch the non-local kernel if there is no work to do.
590        Doing the same for the local kernel is more complicated, since the
591        local part of the force array also depends on the non-local kernel.
592        So to avoid complicating the code and to reduce the risk of bugs,
593        we always call the local kernel and later (not in
594        this function) the stream wait, local f copyback and the f buffer
595        clearing. All these operations, except for the local interaction kernel,
596        are needed for the non-local interactions. The skip of the local kernel
597        call is taken care of later in this function. */
598     if (canSkipNonbondedWork(*nb, iloc))
599     {
600         plist->haveFreshList = false;
601
602         return;
603     }
604
605     if (nbp->useDynamicPruning && plist->haveFreshList)
606     {
607         /* Prunes for rlistOuter and rlistInner, sets plist->haveFreshList=false
608            (that's the way the timing accounting can distinguish between
609            separate prune kernel and combined force+prune).
610          */
611         Nbnxm::gpu_launch_kernel_pruneonly(nb, iloc, 1);
612     }
613
614     if (plist->nsci == 0)
615     {
616         /* Don't launch an empty local kernel (is not allowed with OpenCL).
617          */
618         return;
619     }
620
621     /* beginning of timed nonbonded calculation section */
622     if (bDoTime)
623     {
624         t->interaction[iloc].nb_k.openTimingRegion(stream);
625     }
626
627     /* kernel launch config */
628
629     KernelLaunchConfig config;
630     config.sharedMemorySize = calc_shmem_required_nonbonded(nbp->vdwtype, nb->bPrefetchLjParam);
631     config.stream           = stream;
632     config.blockSize[0]     = c_clSize;
633     config.blockSize[1]     = c_clSize;
634     config.gridSize[0]      = plist->nsci;
635
636     validate_global_work_size(config, 3, nb->dev_info);
637
638     if (debug)
639     {
640         fprintf(debug,
641                 "Non-bonded GPU launch configuration:\n\tLocal work size: %zux%zux%zu\n\t"
642                 "Global work size : %zux%zu\n\t#Super-clusters/clusters: %d/%d (%d)\n",
643                 config.blockSize[0], config.blockSize[1], config.blockSize[2],
644                 config.blockSize[0] * config.gridSize[0], config.blockSize[1] * config.gridSize[1],
645                 plist->nsci * c_numClPerSupercl, c_numClPerSupercl, plist->na_c);
646     }
647
648     fillin_ocl_structures(nbp, &nbparams_params);
649
650     auto*          timingEvent  = bDoTime ? t->interaction[iloc].nb_k.fetchNextEvent() : nullptr;
651     constexpr char kernelName[] = "k_calc_nb";
652     const auto     kernel =
653             select_nbnxn_kernel(nb, nbp->eeltype, nbp->vdwtype, stepWork.computeEnergy,
654                                 (plist->haveFreshList && !nb->timers->interaction[iloc].didPrune));
655
656
657     // The OpenCL kernel takes int as second to last argument because bool is
658     // not supported as a kernel argument type (sizeof(bool) is implementation defined).
659     const int computeFshift = static_cast<int>(stepWork.computeVirial);
660     if (useLjCombRule(nb->nbparam->vdwtype))
661     {
662         const auto kernelArgs = prepareGpuKernelArguments(
663                 kernel, config, &nbparams_params, &adat->xq, &adat->f, &adat->e_lj, &adat->e_el,
664                 &adat->fshift, &adat->lj_comb, &adat->shift_vec, &nbp->nbfp_climg2d, &nbp->nbfp_comb_climg2d,
665                 &nbp->coulomb_tab_climg2d, &plist->sci, &plist->cj4, &plist->excl, &computeFshift);
666
667         launchGpuKernel(kernel, config, timingEvent, kernelName, kernelArgs);
668     }
669     else
670     {
671         const auto kernelArgs = prepareGpuKernelArguments(
672                 kernel, config, &adat->ntypes, &nbparams_params, &adat->xq, &adat->f, &adat->e_lj,
673                 &adat->e_el, &adat->fshift, &adat->atom_types, &adat->shift_vec, &nbp->nbfp_climg2d,
674                 &nbp->nbfp_comb_climg2d, &nbp->coulomb_tab_climg2d, &plist->sci, &plist->cj4,
675                 &plist->excl, &computeFshift);
676         launchGpuKernel(kernel, config, timingEvent, kernelName, kernelArgs);
677     }
678
679     if (bDoTime)
680     {
681         t->interaction[iloc].nb_k.closeTimingRegion(stream);
682     }
683 }
684
685
686 /*! \brief Calculates the amount of shared memory required by the prune kernel.
687  *
688  *  Note that for the sake of simplicity we use the CUDA terminology "shared memory"
689  *  for OpenCL local memory.
690  *
691  * \param[in] num_threads_z cj4 concurrency equal to the number of threads/work items in the 3-rd
692  * dimension. \returns   the amount of local memory in bytes required by the pruning kernel
693  */
694 static inline int calc_shmem_required_prune(const int num_threads_z)
695 {
696     int shmem;
697
698     /* i-atom x in shared memory (for convenience we load all 4 components including q) */
699     shmem = c_numClPerSupercl * c_clSize * sizeof(float) * 4;
700     /* cj in shared memory, for each warp separately
701      * Note: only need to load once per wavefront, but to keep the code simple,
702      * for now we load twice on AMD.
703      */
704     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
705     /* Warp vote, requires one uint per warp/32 threads per block. */
706     shmem += sizeof(cl_uint) * 2 * num_threads_z;
707
708     return shmem;
709 }
710
711 /*! \brief
712  * Launch the pairlist prune only kernel for the given locality.
713  * \p numParts tells in how many parts, i.e. calls the list will be pruned.
714  */
715 void gpu_launch_kernel_pruneonly(gmx_nbnxn_gpu_t* nb, const InteractionLocality iloc, const int numParts)
716 {
717     cl_atomdata_t*   adat    = nb->atdat;
718     cl_nbparam_t*    nbp     = nb->nbparam;
719     cl_plist_t*      plist   = nb->plist[iloc];
720     cl_timers_t*     t       = nb->timers;
721     cl_command_queue stream  = nb->stream[iloc];
722     bool             bDoTime = nb->bDoTime == CL_TRUE;
723
724     if (plist->haveFreshList)
725     {
726         GMX_ASSERT(numParts == 1, "With first pruning we expect 1 part");
727
728         /* Set rollingPruningNumParts to signal that it is not set */
729         plist->rollingPruningNumParts = 0;
730         plist->rollingPruningPart     = 0;
731     }
732     else
733     {
734         if (plist->rollingPruningNumParts == 0)
735         {
736             plist->rollingPruningNumParts = numParts;
737         }
738         else
739         {
740             GMX_ASSERT(numParts == plist->rollingPruningNumParts,
741                        "It is not allowed to change numParts in between list generation steps");
742         }
743     }
744
745     /* Use a local variable for part and update in plist, so we can return here
746      * without duplicating the part increment code.
747      */
748     int part = plist->rollingPruningPart;
749
750     plist->rollingPruningPart++;
751     if (plist->rollingPruningPart >= plist->rollingPruningNumParts)
752     {
753         plist->rollingPruningPart = 0;
754     }
755
756     /* Compute the number of list entries to prune in this pass */
757     int numSciInPart = (plist->nsci - part) / numParts;
758
759     /* Don't launch the kernel if there is no work to do. */
760     if (numSciInPart <= 0)
761     {
762         plist->haveFreshList = false;
763
764         return;
765     }
766
767     GpuRegionTimer* timer = nullptr;
768     if (bDoTime)
769     {
770         timer = &(plist->haveFreshList ? t->interaction[iloc].prune_k : t->interaction[iloc].rollingPrune_k);
771     }
772
773     /* beginning of timed prune calculation section */
774     if (bDoTime)
775     {
776         timer->openTimingRegion(stream);
777     }
778
779     /* Kernel launch config:
780      * - The thread block dimensions match the size of i-clusters, j-clusters,
781      *   and j-cluster concurrency, in x, y, and z, respectively.
782      * - The 1D block-grid contains as many blocks as super-clusters.
783      */
784     int num_threads_z = getOclPruneKernelJ4Concurrency(nb->dev_info->vendor_e);
785
786     /* kernel launch config */
787     KernelLaunchConfig config;
788     config.sharedMemorySize = calc_shmem_required_prune(num_threads_z);
789     config.stream           = stream;
790     config.blockSize[0]     = c_clSize;
791     config.blockSize[1]     = c_clSize;
792     config.blockSize[2]     = num_threads_z;
793     config.gridSize[0]      = numSciInPart;
794
795     validate_global_work_size(config, 3, nb->dev_info);
796
797     if (debug)
798     {
799         fprintf(debug,
800                 "Pruning GPU kernel launch configuration:\n\tLocal work size: %zux%zux%zu\n\t"
801                 "\tGlobal work size: %zux%zu\n\t#Super-clusters/clusters: %d/%d (%d)\n"
802                 "\tShMem: %zu\n",
803                 config.blockSize[0], config.blockSize[1], config.blockSize[2],
804                 config.blockSize[0] * config.gridSize[0], config.blockSize[1] * config.gridSize[1],
805                 plist->nsci * c_numClPerSupercl, c_numClPerSupercl, plist->na_c, config.sharedMemorySize);
806     }
807
808     cl_nbparam_params_t nbparams_params;
809     fillin_ocl_structures(nbp, &nbparams_params);
810
811     auto*          timingEvent  = bDoTime ? timer->fetchNextEvent() : nullptr;
812     constexpr char kernelName[] = "k_pruneonly";
813     const auto     pruneKernel  = selectPruneKernel(nb->kernel_pruneonly, plist->haveFreshList);
814     const auto     kernelArgs   = prepareGpuKernelArguments(pruneKernel, config, &nbparams_params,
815                                                       &adat->xq, &adat->shift_vec, &plist->sci,
816                                                       &plist->cj4, &plist->imask, &numParts, &part);
817     launchGpuKernel(pruneKernel, config, timingEvent, kernelName, kernelArgs);
818
819     if (plist->haveFreshList)
820     {
821         plist->haveFreshList = false;
822         /* Mark that pruning has been done */
823         nb->timers->interaction[iloc].didPrune = true;
824     }
825     else
826     {
827         /* Mark that rolling pruning has been done */
828         nb->timers->interaction[iloc].didRollingPrune = true;
829     }
830
831     if (bDoTime)
832     {
833         timer->closeTimingRegion(stream);
834     }
835 }
836
837 /*! \brief
838  * Launch asynchronously the download of nonbonded forces from the GPU
839  * (and energies/shift forces if required).
840  */
841 void gpu_launch_cpyback(gmx_nbnxn_ocl_t*         nb,
842                         struct nbnxn_atomdata_t* nbatom,
843                         const gmx::StepWorkload& stepWork,
844                         const AtomLocality       aloc)
845 {
846     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
847
848     cl_int gmx_unused cl_error;
849     int               adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
850
851     /* determine interaction locality from atom locality */
852     const InteractionLocality iloc = gpuAtomToInteractionLocality(aloc);
853
854     cl_atomdata_t*   adat    = nb->atdat;
855     cl_timers_t*     t       = nb->timers;
856     bool             bDoTime = nb->bDoTime == CL_TRUE;
857     cl_command_queue stream  = nb->stream[iloc];
858
859     /* don't launch non-local copy-back if there was no non-local work to do */
860     if ((iloc == InteractionLocality::NonLocal) && !haveGpuShortRangeWork(*nb, iloc))
861     {
862         /* TODO An alternative way to signal that non-local work is
863            complete is to use a clEnqueueMarker+clEnqueueBarrier
864            pair. However, the use of bNonLocalStreamActive has the
865            advantage of being local to the host, so probably minimizes
866            overhead. Curiously, for NVIDIA OpenCL with an empty-domain
867            test case, overall simulation performance was higher with
868            the API calls, but this has not been tested on AMD OpenCL,
869            so could be worth considering in future. */
870         nb->bNonLocalStreamActive = CL_FALSE;
871         return;
872     }
873
874     getGpuAtomRange(adat, aloc, &adat_begin, &adat_len);
875
876     /* beginning of timed D2H section */
877     if (bDoTime)
878     {
879         t->xf[aloc].nb_d2h.openTimingRegion(stream);
880     }
881
882     /* With DD the local D2H transfer can only start after the non-local
883        has been launched. */
884     if (iloc == InteractionLocality::Local && nb->bNonLocalStreamActive)
885     {
886         sync_ocl_event(stream, &(nb->nonlocal_done));
887     }
888
889     /* DtoH f */
890     ocl_copy_D2H_async(nbatom->out[0].f.data() + adat_begin * 3, adat->f,
891                        adat_begin * 3 * sizeof(float), (adat_len)*adat->f_elem_size, stream,
892                        bDoTime ? t->xf[aloc].nb_d2h.fetchNextEvent() : nullptr);
893
894     /* kick off work */
895     cl_error = clFlush(stream);
896     assert(CL_SUCCESS == cl_error);
897
898     /* After the non-local D2H is launched the nonlocal_done event can be
899        recorded which signals that the local D2H can proceed. This event is not
900        placed after the non-local kernel because we first need the non-local
901        data back first. */
902     if (iloc == InteractionLocality::NonLocal)
903     {
904         cl_error = clEnqueueMarkerWithWaitList(stream, 0, nullptr, &(nb->nonlocal_done));
905         assert(CL_SUCCESS == cl_error);
906         nb->bNonLocalStreamActive = CL_TRUE;
907     }
908
909     /* only transfer energies in the local stream */
910     if (iloc == InteractionLocality::Local)
911     {
912         /* DtoH fshift when virial is needed */
913         if (stepWork.computeVirial)
914         {
915             ocl_copy_D2H_async(nb->nbst.fshift, adat->fshift, 0, SHIFTS * adat->fshift_elem_size,
916                                stream, bDoTime ? t->xf[aloc].nb_d2h.fetchNextEvent() : nullptr);
917         }
918
919         /* DtoH energies */
920         if (stepWork.computeEnergy)
921         {
922             ocl_copy_D2H_async(nb->nbst.e_lj, adat->e_lj, 0, sizeof(float), stream,
923                                bDoTime ? t->xf[aloc].nb_d2h.fetchNextEvent() : nullptr);
924
925             ocl_copy_D2H_async(nb->nbst.e_el, adat->e_el, 0, sizeof(float), stream,
926                                bDoTime ? t->xf[aloc].nb_d2h.fetchNextEvent() : nullptr);
927         }
928     }
929
930     if (bDoTime)
931     {
932         t->xf[aloc].nb_d2h.closeTimingRegion(stream);
933     }
934 }
935
936
937 /*! \brief Selects the Ewald kernel type, analytical or tabulated, single or twin cut-off. */
938 int nbnxn_gpu_pick_ewald_kernel_type(const interaction_const_t& ic)
939 {
940     bool bTwinCut = (ic.rcoulomb != ic.rvdw);
941     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
942     int  kernel_type;
943
944     /* Benchmarking/development environment variables to force the use of
945        analytical or tabulated Ewald kernel. */
946     bForceAnalyticalEwald = (getenv("GMX_OCL_NB_ANA_EWALD") != nullptr);
947     bForceTabulatedEwald  = (getenv("GMX_OCL_NB_TAB_EWALD") != nullptr);
948
949     if (bForceAnalyticalEwald && bForceTabulatedEwald)
950     {
951         gmx_incons(
952                 "Both analytical and tabulated Ewald OpenCL non-bonded kernels "
953                 "requested through environment variables.");
954     }
955
956     /* OpenCL: By default, use analytical Ewald
957      * TODO: tabulated does not work, it needs fixing, see init_nbparam() in nbnxn_ocl_data_mgmt.cpp
958      *
959      * TODO: decide if dev_info parameter should be added to recognize NVIDIA CC>=3.0 devices.
960      *
961      */
962     /* By default use analytical Ewald. */
963     bUseAnalyticalEwald = true;
964     if (bForceAnalyticalEwald)
965     {
966         if (debug)
967         {
968             fprintf(debug, "Using analytical Ewald OpenCL kernels\n");
969         }
970     }
971     else if (bForceTabulatedEwald)
972     {
973         bUseAnalyticalEwald = false;
974
975         if (debug)
976         {
977             fprintf(debug, "Using tabulated Ewald OpenCL kernels\n");
978         }
979     }
980
981     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
982        forces it (use it for debugging/benchmarking only). */
983     if (!bTwinCut && (getenv("GMX_OCL_NB_EWALD_TWINCUT") == nullptr))
984     {
985         kernel_type = bUseAnalyticalEwald ? eelOclEWALD_ANA : eelOclEWALD_TAB;
986     }
987     else
988     {
989         kernel_type = bUseAnalyticalEwald ? eelOclEWALD_ANA_TWIN : eelOclEWALD_TAB_TWIN;
990     }
991
992     return kernel_type;
993 }
994
995 } // namespace Nbnxm