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