9a1a40ff0baaca549a6e9ddc46b1a32b81b4b7b6
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_ocl / nbnxn_ocl.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017, 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 nbnxn_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_mdlib
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/oclutils.h"
73 #include "gromacs/hardware/hw_info.h"
74 #include "gromacs/mdlib/force_flags.h"
75 #include "gromacs/mdlib/nb_verlet.h"
76 #include "gromacs/mdlib/nbnxn_consts.h"
77 #include "gromacs/mdlib/nbnxn_gpu.h"
78 #include "gromacs/mdlib/nbnxn_gpu_common.h"
79 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
80 #include "gromacs/mdlib/nbnxn_pairlist.h"
81 #include "gromacs/pbcutil/ishift.h"
82 #include "gromacs/timing/gpu_timing.h"
83 #include "gromacs/utility/cstringutil.h"
84 #include "gromacs/utility/fatalerror.h"
85 #include "gromacs/utility/gmxassert.h"
86
87 #include "nbnxn_ocl_internal.h"
88 #include "nbnxn_ocl_types.h"
89
90 #if defined TEXOBJ_SUPPORTED && __CUDA_ARCH__ >= 300
91 #define USE_TEXOBJ
92 #endif
93
94 /*! \brief Convenience constants */
95 //@{
96 static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
97 static const int c_clSize          = c_nbnxnGpuClusterSize;
98 //@}
99
100 /*! \brief Always/never run the energy/pruning kernels -- only for benchmarking purposes */
101 //@{
102 static bool always_ener  = (getenv("GMX_GPU_ALWAYS_ENER") != NULL);
103 static bool never_ener   = (getenv("GMX_GPU_NEVER_ENER") != NULL);
104 static bool always_prune = (getenv("GMX_GPU_ALWAYS_PRUNE") != NULL);
105 //@}
106
107 /* Uncomment this define to enable kernel debugging */
108 //#define DEBUG_OCL
109
110 /*! \brief Specifies which kernel run to debug */
111 #define DEBUG_RUN_STEP 2
112
113 /*! \brief Validates the input global work size parameter.
114  */
115 static inline void validate_global_work_size(size_t *global_work_size, int work_dim, const gmx_device_info_t *dinfo)
116 {
117     cl_uint device_size_t_size_bits;
118     cl_uint host_size_t_size_bits;
119
120     assert(dinfo);
121
122     /* Each component of a global_work_size must not exceed the range given by the
123        sizeof(device size_t) for the device on which the kernel execution will
124        be enqueued. See:
125        https://www.khronos.org/registry/cl/sdk/1.0/docs/man/xhtml/clEnqueueNDRangeKernel.html
126      */
127     device_size_t_size_bits = dinfo->adress_bits;
128     host_size_t_size_bits   = (cl_uint)(sizeof(size_t) * 8);
129
130     /* If sizeof(host size_t) <= sizeof(device size_t)
131             => global_work_size components will always be valid
132        else
133             => get device limit for global work size and
134             compare it against each component of global_work_size.
135      */
136     if (host_size_t_size_bits > device_size_t_size_bits)
137     {
138         size_t device_limit;
139
140         device_limit = (((size_t)1) << device_size_t_size_bits) - 1;
141
142         for (int i = 0; i < work_dim; i++)
143         {
144             if (global_work_size[i] > device_limit)
145             {
146                 gmx_fatal(FARGS, "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 (%d > %d)!",
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 nbnxn_cuda_types.h.
161  */
162
163 /*! \brief Force-only kernel function names. */
164 static const char* nb_kfunc_noener_noprune_ptr[eelOclNR][evdwOclNR] =
165 {
166     { "nbnxn_kernel_ElecCut_VdwLJ_F_opencl",            "nbnxn_kernel_ElecCut_VdwLJCombGeom_F_opencl",            "nbnxn_kernel_ElecCut_VdwLJCombLB_F_opencl",            "nbnxn_kernel_ElecCut_VdwLJFsw_F_opencl",            "nbnxn_kernel_ElecCut_VdwLJPsw_F_opencl",            "nbnxn_kernel_ElecCut_VdwLJEwCombGeom_F_opencl",            "nbnxn_kernel_ElecCut_VdwLJEwCombLB_F_opencl"            },
167     { "nbnxn_kernel_ElecRF_VdwLJ_F_opencl",             "nbnxn_kernel_ElecRF_VdwLJCombGeom_F_opencl",             "nbnxn_kernel_ElecRF_VdwLJCombLB_F_opencl",             "nbnxn_kernel_ElecRF_VdwLJFsw_F_opencl",             "nbnxn_kernel_ElecRF_VdwLJPsw_F_opencl",             "nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_opencl",             "nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_opencl"             },
168     { "nbnxn_kernel_ElecEwQSTab_VdwLJ_F_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_F_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_F_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_F_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_F_opencl"        },
169     { "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_F_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_F_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_F_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_F_opencl" },
170     { "nbnxn_kernel_ElecEw_VdwLJ_F_opencl",             "nbnxn_kernel_ElecEw_VdwLJCombGeom_F_opencl",             "nbnxn_kernel_ElecEw_VdwLJCombLB_F_opencl",             "nbnxn_kernel_ElecEw_VdwLJFsw_F_opencl",             "nbnxn_kernel_ElecEw_VdwLJPsw_F_opencl",             "nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_opencl",             "nbnxn_kernel_ElecEw_VdwLJEwCombLB_F_opencl"             },
171     { "nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_F_opencl"      }
172 };
173
174 /*! \brief Force + energy kernel function pointers. */
175 static const char* nb_kfunc_ener_noprune_ptr[eelOclNR][evdwOclNR] =
176 {
177     { "nbnxn_kernel_ElecCut_VdwLJ_VF_opencl",            "nbnxn_kernel_ElecCut_VdwLJCombGeom_VF_opencl",            "nbnxn_kernel_ElecCut_VdwLJCombLB_VF_opencl",            "nbnxn_kernel_ElecCut_VdwLJFsw_VF_opencl",            "nbnxn_kernel_ElecCut_VdwLJPsw_VF_opencl",            "nbnxn_kernel_ElecCut_VdwLJEwCombGeom_VF_opencl",            "nbnxn_kernel_ElecCut_VdwLJEwCombLB_VF_opencl"            },
178     { "nbnxn_kernel_ElecRF_VdwLJ_VF_opencl",             "nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_opencl",             "nbnxn_kernel_ElecRF_VdwLJCombLB_VF_opencl",             "nbnxn_kernel_ElecRF_VdwLJFsw_VF_opencl",             "nbnxn_kernel_ElecRF_VdwLJPsw_VF_opencl",             "nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_opencl",             "nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_opencl"             },
179     { "nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_VF_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_VF_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_VF_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_VF_opencl"        },
180     { "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_VF_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_VF_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_VF_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_VF_opencl" },
181     { "nbnxn_kernel_ElecEw_VdwLJ_VF_opencl",             "nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_opencl",             "nbnxn_kernel_ElecEw_VdwLJCombLB_VF_opencl",             "nbnxn_kernel_ElecEw_VdwLJFsw_VF_opencl",             "nbnxn_kernel_ElecEw_VdwLJPsw_VF_opencl",             "nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_opencl",             "nbnxn_kernel_ElecEw_VdwLJEwCombLB_VF_opencl"             },
182     { "nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_VF_opencl"      }
183 };
184
185 /*! \brief Force + pruning kernel function pointers. */
186 static const char* nb_kfunc_noener_prune_ptr[eelOclNR][evdwOclNR] =
187 {
188     { "nbnxn_kernel_ElecCut_VdwLJ_F_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJCombGeom_F_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJCombLB_F_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJFsw_F_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJPsw_F_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJEwCombGeom_F_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJEwCombLB_F_prune_opencl"             },
189     { "nbnxn_kernel_ElecRF_VdwLJ_F_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJCombGeom_F_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJCombLB_F_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJFsw_F_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJPsw_F_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_prune_opencl"              },
190     { "nbnxn_kernel_ElecEwQSTab_VdwLJ_F_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_F_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_F_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_F_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_F_prune_opencl"         },
191     { "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_F_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_F_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_F_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_F_prune_opencl"  },
192     { "nbnxn_kernel_ElecEw_VdwLJ_F_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJCombGeom_F_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJCombLB_F_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJFsw_F_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJPsw_F_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJEwCombLB_F_prune_opencl"              },
193     { "nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_F_prune_opencl"       }
194 };
195
196 /*! \brief Force + energy + pruning kernel function pointers. */
197 static const char* nb_kfunc_ener_prune_ptr[eelOclNR][evdwOclNR] =
198 {
199     { "nbnxn_kernel_ElecCut_VdwLJ_VF_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJCombGeom_VF_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJCombLB_VF_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJFsw_VF_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJPsw_VF_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJEwCombGeom_VF_prune_opencl",            "nbnxn_kernel_ElecCut_VdwLJEwCombLB_VF_prune_opencl"            },
200     { "nbnxn_kernel_ElecRF_VdwLJ_VF_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJCombLB_VF_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJFsw_VF_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJPsw_VF_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_prune_opencl",             "nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_prune_opencl"             },
201     { "nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_VF_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_VF_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_VF_prune_opencl",        "nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_VF_prune_opencl"        },
202     { "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_VF_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_VF_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_VF_prune_opencl", "nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_VF_prune_opencl" },
203     { "nbnxn_kernel_ElecEw_VdwLJ_VF_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJCombLB_VF_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJFsw_VF_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJPsw_VF_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_prune_opencl",             "nbnxn_kernel_ElecEw_VdwLJEwCombLB_VF_prune_opencl"             },
204     { "nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_prune_opencl",      "nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_VF_prune_opencl"      }
205 };
206
207 /*! \brief Return a pointer to the prune kernel version to be executed at the current invocation.
208  *
209  * \param[in] kernel_pruneonly  array of prune kernel objects
210  * \param[in] firstPrunePass    true if the first pruning pass is being executed
211  */
212 static inline cl_kernel selectPruneKernel(cl_kernel kernel_pruneonly[],
213                                           bool      firstPrunePass)
214 {
215     cl_kernel  *kernelPtr;
216
217     if (firstPrunePass)
218     {
219         kernelPtr = &(kernel_pruneonly[epruneFirst]);
220     }
221     else
222     {
223         kernelPtr = &(kernel_pruneonly[epruneRolling]);
224     }
225     // TODO: consider creating the prune kernel object here to avoid a
226     // clCreateKernel for the rolling prune kernel if this is not needed.
227     return *kernelPtr;
228 }
229
230 /*! \brief Return a pointer to the kernel version to be executed at the current step.
231  *  OpenCL kernel objects are cached in nb. If the requested kernel is not
232  *  found in the cache, it will be created and the cache will be updated.
233  */
234 static inline cl_kernel select_nbnxn_kernel(gmx_nbnxn_ocl_t   *nb,
235                                             int                eeltype,
236                                             int                evdwtype,
237                                             bool               bDoEne,
238                                             bool               bDoPrune)
239 {
240     const char* kernel_name_to_run;
241     cl_kernel  *kernel_ptr;
242     cl_int      cl_error;
243
244     assert(eeltype  < eelOclNR);
245     assert(evdwtype < evdwOclNR);
246
247     if (bDoEne)
248     {
249         if (bDoPrune)
250         {
251             kernel_name_to_run = nb_kfunc_ener_prune_ptr[eeltype][evdwtype];
252             kernel_ptr         = &(nb->kernel_ener_prune_ptr[eeltype][evdwtype]);
253         }
254         else
255         {
256             kernel_name_to_run = nb_kfunc_ener_noprune_ptr[eeltype][evdwtype];
257             kernel_ptr         = &(nb->kernel_ener_noprune_ptr[eeltype][evdwtype]);
258         }
259     }
260     else
261     {
262         if (bDoPrune)
263         {
264             kernel_name_to_run = nb_kfunc_noener_prune_ptr[eeltype][evdwtype];
265             kernel_ptr         = &(nb->kernel_noener_prune_ptr[eeltype][evdwtype]);
266         }
267         else
268         {
269             kernel_name_to_run = nb_kfunc_noener_noprune_ptr[eeltype][evdwtype];
270             kernel_ptr         = &(nb->kernel_noener_noprune_ptr[eeltype][evdwtype]);
271         }
272     }
273
274     if (NULL == kernel_ptr[0])
275     {
276         *kernel_ptr = clCreateKernel(nb->dev_rundata->program, kernel_name_to_run, &cl_error);
277         assert(cl_error == CL_SUCCESS);
278     }
279     // TODO: handle errors
280
281     return *kernel_ptr;
282 }
283
284 /*! \brief Calculates the amount of shared memory required by the nonbonded kernel in use.
285  */
286 static inline int calc_shmem_required_nonbonded(int  vdwType,
287                                                 bool bPrefetchLjParam)
288 {
289     int shmem;
290
291     /* size of shmem (force-buffers/xq/atom type preloading) */
292     /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
293     /* i-atom x+q in shared memory */
294     shmem  = c_numClPerSupercl * c_clSize * sizeof(float) * 4; /* xqib */
295     /* cj in shared memory, for both warps separately */
296     shmem += 2 * c_nbnxnGpuJgroupSize * sizeof(int);           /* cjs  */
297     if (bPrefetchLjParam)
298     {
299         if (useLjCombRule(vdwType))
300         {
301             /* i-atom LJ combination parameters in shared memory */
302             shmem += c_numClPerSupercl * c_clSize * 2*sizeof(float); /* atib abused for ljcp, float2 */
303         }
304         else
305         {
306             /* i-atom types in shared memory */
307             shmem += c_numClPerSupercl * c_clSize * sizeof(int); /* atib */
308         }
309     }
310     /* force reduction buffers in shared memory */
311     shmem += c_clSize * c_clSize * 3 * sizeof(float);    /* f_buf */
312     /* Warp vote. In fact it must be * number of warps in block.. */
313     shmem += sizeof(cl_uint) * 2;                        /* warp_any */
314     return shmem;
315 }
316
317 /*! \brief Initializes data structures that are going to be sent to the OpenCL device.
318  *
319  *  The device can't use the same data structures as the host for two main reasons:
320  *  - OpenCL restrictions (pointers are not accepted inside data structures)
321  *  - some host side fields are not needed for the OpenCL kernels.
322  *
323  *  This function is called before the launch of both nbnxn and prune kernels.
324  */
325 static void fillin_ocl_structures(cl_nbparam_t        *nbp,
326                                   cl_nbparam_params_t *nbparams_params)
327 {
328     nbparams_params->coulomb_tab_scale = nbp->coulomb_tab_scale;
329     nbparams_params->c_rf              = nbp->c_rf;
330     nbparams_params->dispersion_shift  = nbp->dispersion_shift;
331     nbparams_params->eeltype           = nbp->eeltype;
332     nbparams_params->epsfac            = nbp->epsfac;
333     nbparams_params->ewaldcoeff_lj     = nbp->ewaldcoeff_lj;
334     nbparams_params->ewald_beta        = nbp->ewald_beta;
335     nbparams_params->rcoulomb_sq       = nbp->rcoulomb_sq;
336     nbparams_params->repulsion_shift   = nbp->repulsion_shift;
337     nbparams_params->rlistOuter_sq     = nbp->rlistOuter_sq;
338     nbparams_params->rvdw_sq           = nbp->rvdw_sq;
339     nbparams_params->rlistInner_sq     = nbp->rlistInner_sq;
340     nbparams_params->rvdw_switch       = nbp->rvdw_switch;
341     nbparams_params->sh_ewald          = nbp->sh_ewald;
342     nbparams_params->sh_lj_ewald       = nbp->sh_lj_ewald;
343     nbparams_params->two_k_rf          = nbp->two_k_rf;
344     nbparams_params->vdwtype           = nbp->vdwtype;
345     nbparams_params->vdw_switch        = nbp->vdw_switch;
346 }
347
348 /*! \brief Enqueues a wait for event completion.
349  *
350  * Then it releases the event and sets it to 0.
351  * Don't use this function when more than one wait will be issued for the event.
352  * Equivalent to Cuda Stream Sync. */
353 static void sync_ocl_event(cl_command_queue stream, cl_event *ocl_event)
354 {
355     cl_int gmx_unused cl_error;
356
357     /* Enqueue wait */
358 #ifdef CL_VERSION_1_2
359     cl_error = clEnqueueBarrierWithWaitList(stream, 1, ocl_event, NULL);
360 #else
361     cl_error = clEnqueueWaitForEvents(stream, 1, ocl_event);
362 #endif
363
364     GMX_RELEASE_ASSERT(CL_SUCCESS == cl_error, ocl_get_error_string(cl_error).c_str());
365
366     /* Release event and reset it to 0. It is ok to release it as enqueuewaitforevents performs implicit retain for events. */
367     cl_error = clReleaseEvent(*ocl_event);
368     assert(CL_SUCCESS == cl_error);
369     *ocl_event = 0;
370 }
371
372 /*! \brief Launch GPU kernel
373
374    As we execute nonbonded workload in separate queues, before launching
375    the kernel we need to make sure that he following operations have completed:
376    - atomdata allocation and related H2D transfers (every nstlist step);
377    - pair list H2D transfer (every nstlist step);
378    - shift vector H2D transfer (every nstlist step);
379    - force (+shift force and energy) output clearing (every step).
380
381    These operations are issued in the local queue at the beginning of the step
382    and therefore always complete before the local kernel launch. The non-local
383    kernel is launched after the local on the same device/context, so this is
384    inherently scheduled after the operations in the local stream (including the
385    above "misc_ops").
386    However, for the sake of having a future-proof implementation, we use the
387    misc_ops_done event to record the point in time when the above  operations
388    are finished and synchronize with this event in the non-local stream.
389  */
390 void nbnxn_gpu_launch_kernel(gmx_nbnxn_ocl_t               *nb,
391                              const struct nbnxn_atomdata_t *nbatom,
392                              int                            flags,
393                              int                            iloc)
394 {
395     cl_int               cl_error;
396     int                  adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
397     /* OpenCL kernel launch-related stuff */
398     int                  shmem;
399     size_t               local_work_size[3], global_work_size[3];
400     cl_kernel            nb_kernel = NULL; /* fn pointer to the nonbonded kernel */
401
402     cl_atomdata_t       *adat    = nb->atdat;
403     cl_nbparam_t        *nbp     = nb->nbparam;
404     cl_plist_t          *plist   = nb->plist[iloc];
405     cl_timers_t         *t       = nb->timers;
406     cl_command_queue     stream  = nb->stream[iloc];
407
408     bool                 bCalcEner   = flags & GMX_FORCE_ENERGY;
409     int                  bCalcFshift = flags & GMX_FORCE_VIRIAL;
410     bool                 bDoTime     = nb->bDoTime;
411     cl_uint              arg_no;
412
413     cl_nbparam_params_t  nbparams_params;
414 #ifdef DEBUG_OCL
415     float              * debug_buffer_h;
416     size_t               debug_buffer_size;
417 #endif
418
419     /* turn energy calculation always on/off (for debugging/testing only) */
420     bCalcEner = (bCalcEner || always_ener) && !never_ener;
421
422     /* Don't launch the non-local kernel if there is no work to do.
423        Doing the same for the local kernel is more complicated, since the
424        local part of the force array also depends on the non-local kernel.
425        So to avoid complicating the code and to reduce the risk of bugs,
426        we always call the local kernel, the local x+q copy and later (not in
427        this function) the stream wait, local f copyback and the f buffer
428        clearing. All these operations, except for the local interaction kernel,
429        are needed for the non-local interactions. The skip of the local kernel
430        call is taken care of later in this function. */
431     if (canSkipWork(nb, iloc))
432     {
433         plist->haveFreshList = false;
434
435         return;
436     }
437
438     /* calculate the atom data index range based on locality */
439     if (LOCAL_I(iloc))
440     {
441         adat_begin  = 0;
442         adat_len    = adat->natoms_local;
443     }
444     else
445     {
446         adat_begin  = adat->natoms_local;
447         adat_len    = adat->natoms - adat->natoms_local;
448     }
449
450     /* beginning of timed HtoD section */
451     if (bDoTime)
452     {
453         t->nb_h2d[iloc].openTimingRegion(stream);
454     }
455
456     /* HtoD x, q */
457     ocl_copy_H2D_async(adat->xq, nbatom->x + adat_begin * 4, adat_begin*sizeof(float)*4,
458                        adat_len * sizeof(float) * 4, stream, bDoTime ? t->nb_h2d[iloc].fetchNextEvent() : nullptr);
459
460     if (bDoTime)
461     {
462         t->nb_h2d[iloc].closeTimingRegion(stream);
463     }
464
465     /* When we get here all misc operations issues in the local stream as well as
466        the local xq H2D are done,
467        so we record that in the local stream and wait for it in the nonlocal one. */
468     if (nb->bUseTwoStreams)
469     {
470         if (iloc == eintLocal)
471         {
472 #ifdef CL_VERSION_1_2
473             cl_error = clEnqueueMarkerWithWaitList(stream, 0, NULL, &(nb->misc_ops_and_local_H2D_done));
474 #else
475             cl_error = clEnqueueMarker(stream, &(nb->misc_ops_and_local_H2D_done));
476 #endif
477             assert(CL_SUCCESS == cl_error);
478
479             /* Based on the v1.2 section 5.13 of the OpenCL spec, a flush is needed
480              * in the local stream in order to be able to sync with the above event
481              * from the non-local stream.
482              */
483             cl_error = clFlush(stream);
484             assert(CL_SUCCESS == cl_error);
485         }
486         else
487         {
488             sync_ocl_event(stream, &(nb->misc_ops_and_local_H2D_done));
489         }
490     }
491
492     if (nbp->useDynamicPruning && plist->haveFreshList)
493     {
494         /* Prunes for rlistOuter and rlistInner, sets plist->haveFreshList=false
495            (that's the way the timing accounting can distinguish between
496            separate prune kernel and combined force+prune).
497          */
498         nbnxn_gpu_launch_kernel_pruneonly(nb, iloc, 1);
499     }
500
501     if (plist->nsci == 0)
502     {
503         /* Don't launch an empty local kernel (is not allowed with OpenCL).
504          * TODO: Separate H2D and kernel launch into separate functions.
505          */
506         return;
507     }
508
509     /* beginning of timed nonbonded calculation section */
510     if (bDoTime)
511     {
512         t->nb_k[iloc].openTimingRegion(stream);
513     }
514
515     /* get the pointer to the kernel flavor we need to use */
516     nb_kernel = select_nbnxn_kernel(nb,
517                                     nbp->eeltype,
518                                     nbp->vdwtype,
519                                     bCalcEner,
520                                     (plist->haveFreshList && !nb->timers->didPrune[iloc]) || always_prune);
521
522     /* kernel launch config */
523     local_work_size[0] = c_clSize;
524     local_work_size[1] = c_clSize;
525     local_work_size[2] = 1;
526
527     global_work_size[0] = plist->nsci * local_work_size[0];
528     global_work_size[1] = 1 * local_work_size[1];
529     global_work_size[2] = 1 * local_work_size[2];
530
531     validate_global_work_size(global_work_size, 3, nb->dev_info);
532
533     shmem     = calc_shmem_required_nonbonded(nbp->vdwtype, nb->bPrefetchLjParam);
534
535 #ifdef DEBUG_OCL
536     {
537         static int run_step = 1;
538
539         if (DEBUG_RUN_STEP == run_step)
540         {
541             debug_buffer_size = global_work_size[0] * global_work_size[1] * global_work_size[2] * sizeof(float);
542             debug_buffer_h    = (float*)calloc(1, debug_buffer_size);
543             assert(NULL != debug_buffer_h);
544
545             if (NULL == nb->debug_buffer)
546             {
547                 nb->debug_buffer = clCreateBuffer(nb->dev_rundata->context, CL_MEM_READ_WRITE | CL_MEM_COPY_HOST_PTR,
548                                                   debug_buffer_size, debug_buffer_h, &cl_error);
549
550                 assert(CL_SUCCESS == cl_error);
551             }
552         }
553
554         run_step++;
555     }
556 #endif
557     if (debug)
558     {
559         fprintf(debug, "Non-bonded GPU launch configuration:\n\tLocal work size: %dx%dx%d\n\t"
560                 "Global work size : %dx%d\n\t#Super-clusters/clusters: %d/%d (%d)\n",
561                 (int)(local_work_size[0]), (int)(local_work_size[1]), (int)(local_work_size[2]),
562                 (int)(global_work_size[0]), (int)(global_work_size[1]), plist->nsci*c_numClPerSupercl,
563                 c_numClPerSupercl, plist->na_c);
564     }
565
566     fillin_ocl_structures(nbp, &nbparams_params);
567
568     arg_no    = 0;
569     cl_error  = CL_SUCCESS;
570     if (!useLjCombRule(nb->nbparam->vdwtype))
571     {
572         cl_error  = clSetKernelArg(nb_kernel, arg_no++, sizeof(int), &(adat->ntypes));
573     }
574     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(nbparams_params), &(nbparams_params));
575     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(adat->xq));
576     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(adat->f));
577     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(adat->e_lj));
578     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(adat->e_el));
579     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(adat->fshift));
580     if (useLjCombRule(nb->nbparam->vdwtype))
581     {
582         cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(adat->lj_comb));
583     }
584     else
585     {
586         cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(adat->atom_types));
587     }
588     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(adat->shift_vec));
589     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(nbp->nbfp_climg2d));
590     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(nbp->nbfp_comb_climg2d));
591     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(nbp->coulomb_tab_climg2d));
592     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(plist->sci));
593     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(plist->cj4));
594     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(plist->excl));
595     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(int), &bCalcFshift);
596     cl_error |= clSetKernelArg(nb_kernel, arg_no++, shmem, NULL);
597     cl_error |= clSetKernelArg(nb_kernel, arg_no++, sizeof(cl_mem), &(nb->debug_buffer));
598
599     assert(cl_error == CL_SUCCESS);
600
601     if (cl_error)
602     {
603         printf("OpenCL error: %s\n", ocl_get_error_string(cl_error).c_str());
604     }
605     cl_error = clEnqueueNDRangeKernel(stream, nb_kernel, 3, NULL, global_work_size, local_work_size, 0, NULL, bDoTime ? t->nb_k[iloc].fetchNextEvent() : nullptr);
606     assert(cl_error == CL_SUCCESS);
607
608     if (bDoTime)
609     {
610         t->nb_k[iloc].closeTimingRegion(stream);
611     }
612
613 #ifdef DEBUG_OCL
614     {
615         static int run_step = 1;
616
617         if (DEBUG_RUN_STEP == run_step)
618         {
619             FILE *pf;
620             char  file_name[256] = {0};
621
622             ocl_copy_D2H_async(debug_buffer_h, nb->debug_buffer, 0,
623                                debug_buffer_size, stream, NULL);
624
625             // Make sure all data has been transfered back from device
626             clFinish(stream);
627
628             printf("\nWriting debug_buffer to debug_buffer_ocl.txt...");
629
630             sprintf(file_name, "debug_buffer_ocl_%d.txt", DEBUG_RUN_STEP);
631             pf = fopen(file_name, "wt");
632             assert(pf != NULL);
633
634             fprintf(pf, "%20s", "");
635             for (int j = 0; j < global_work_size[0]; j++)
636             {
637                 char label[20];
638                 sprintf(label, "(wIdx=%2d thIdx=%2d)", j / local_work_size[0], j % local_work_size[0]);
639                 fprintf(pf, "%20s", label);
640             }
641
642             for (int i = 0; i < global_work_size[1]; i++)
643             {
644                 char label[20];
645                 sprintf(label, "(wIdy=%2d thIdy=%2d)", i / local_work_size[1], i % local_work_size[1]);
646                 fprintf(pf, "\n%20s", label);
647
648                 for (int j = 0; j < global_work_size[0]; j++)
649                 {
650                     fprintf(pf, "%20.5f", debug_buffer_h[i * global_work_size[0] + j]);
651                 }
652
653                 //fprintf(pf, "\n");
654             }
655
656             fclose(pf);
657
658             printf(" done.\n");
659
660
661             free(debug_buffer_h);
662             debug_buffer_h = NULL;
663         }
664
665         run_step++;
666     }
667 #endif
668 }
669
670
671 /*! \brief Calculates the amount of shared memory required by the prune kernel.
672  *
673  *  Note that for the sake of simplicity we use the CUDA terminology "shared memory"
674  *  for OpenCL local memory.
675  *
676  * \param[in] num_threads_z cj4 concurrency equal to the number of threads/work items in the 3-rd dimension.
677  * \returns   the amount of local memory in bytes required by the pruning kernel
678  */
679 static inline int calc_shmem_required_prune(const int num_threads_z)
680 {
681     int shmem;
682
683     /* i-atom x in shared memory (for convenience we load all 4 components including q) */
684     shmem  = c_numClPerSupercl * c_clSize * sizeof(float)*4;
685     /* cj in shared memory, for each warp separately */
686     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
687     /* Warp vote, requires one uint per warp/32 threads per block. */
688     shmem += sizeof(cl_uint) * 2*num_threads_z;
689
690     return shmem;
691 }
692
693 void nbnxn_gpu_launch_kernel_pruneonly(gmx_nbnxn_gpu_t       *nb,
694                                        int                    iloc,
695                                        int                    numParts)
696 {
697     cl_int               cl_error;
698
699     cl_atomdata_t       *adat    = nb->atdat;
700     cl_nbparam_t        *nbp     = nb->nbparam;
701     cl_plist_t          *plist   = nb->plist[iloc];
702     cl_timers_t         *t       = nb->timers;
703     cl_command_queue     stream  = nb->stream[iloc];
704     bool                 bDoTime = nb->bDoTime;
705
706     if (plist->haveFreshList)
707     {
708         GMX_ASSERT(numParts == 1, "With first pruning we expect 1 part");
709
710         /* Set rollingPruningNumParts to signal that it is not set */
711         plist->rollingPruningNumParts = 0;
712         plist->rollingPruningPart     = 0;
713     }
714     else
715     {
716         if (plist->rollingPruningNumParts == 0)
717         {
718             plist->rollingPruningNumParts = numParts;
719         }
720         else
721         {
722             GMX_ASSERT(numParts == plist->rollingPruningNumParts, "It is not allowed to change numParts in between list generation steps");
723         }
724     }
725
726     /* Use a local variable for part and update in plist, so we can return here
727      * without duplicating the part increment code.
728      */
729     int part = plist->rollingPruningPart;
730
731     plist->rollingPruningPart++;
732     if (plist->rollingPruningPart >= plist->rollingPruningNumParts)
733     {
734         plist->rollingPruningPart = 0;
735     }
736
737     /* Compute the number of list entries to prune in this pass */
738     int numSciInPart = (plist->nsci - part)/numParts;
739
740     /* Don't launch the kernel if there is no work to do. */
741     if (numSciInPart <= 0)
742     {
743         plist->haveFreshList = false;
744
745         return;
746     }
747
748     GpuRegionTimer *timer = nullptr;
749     if (bDoTime)
750     {
751         timer = &(plist->haveFreshList ? t->prune_k[iloc] : t->rollingPrune_k[iloc]);
752     }
753
754     /* beginning of timed prune calculation section */
755     if (bDoTime)
756     {
757         timer->openTimingRegion(stream);
758     }
759
760     /* Kernel launch config:
761      * - The thread block dimensions match the size of i-clusters, j-clusters,
762      *   and j-cluster concurrency, in x, y, and z, respectively.
763      * - The 1D block-grid contains as many blocks as super-clusters.
764      */
765     int       num_threads_z = getOclPruneKernelJ4Concurrency(nb->dev_info->vendor_e);
766     cl_kernel pruneKernel   = selectPruneKernel(nb->kernel_pruneonly, plist->haveFreshList);
767
768     /* kernel launch config */
769     size_t  local_work_size[3], global_work_size[3];
770     local_work_size[0] = c_clSize;
771     local_work_size[1] = c_clSize;
772     local_work_size[2] = num_threads_z;
773
774     global_work_size[0] = numSciInPart * local_work_size[0];
775     global_work_size[1] = 1 * local_work_size[1];
776     global_work_size[2] = 1 * local_work_size[2];
777
778     validate_global_work_size(global_work_size, 3, nb->dev_info);
779
780     int shmem = calc_shmem_required_prune(num_threads_z);
781
782     if (debug)
783     {
784         fprintf(debug, "Pruning GPU kernel launch configuration:\n\tLocal work size: %dx%dx%d\n\t"
785                 "\tGlobal work size: %dx%d\n\t#Super-clusters/clusters: %d/%d (%d)\n"
786                 "\tShMem: %d\n",
787                 (int)(local_work_size[0]), (int)(local_work_size[1]), (int)(local_work_size[2]),
788                 (int)(global_work_size[0]), (int)(global_work_size[1]), plist->nsci*c_numClPerSupercl,
789                 c_numClPerSupercl, plist->na_c, shmem);
790     }
791
792     cl_nbparam_params_t  nbparams_params;
793     fillin_ocl_structures(nbp, &nbparams_params);
794
795     cl_uint  arg_no = 0;
796     cl_error = CL_SUCCESS;
797
798     cl_error |= clSetKernelArg(pruneKernel, arg_no++, sizeof(nbparams_params), &(nbparams_params));
799     cl_error |= clSetKernelArg(pruneKernel, arg_no++, sizeof(cl_mem), &(adat->xq));
800     cl_error |= clSetKernelArg(pruneKernel, arg_no++, sizeof(cl_mem), &(adat->shift_vec));
801     cl_error |= clSetKernelArg(pruneKernel, arg_no++, sizeof(cl_mem), &(plist->sci));
802     cl_error |= clSetKernelArg(pruneKernel, arg_no++, sizeof(cl_mem), &(plist->cj4));
803     cl_error |= clSetKernelArg(pruneKernel, arg_no++, sizeof(cl_mem), &(plist->imask));
804     cl_error |= clSetKernelArg(pruneKernel, arg_no++, sizeof(int), &(numParts));
805     cl_error |= clSetKernelArg(pruneKernel, arg_no++, sizeof(int), &(part));
806     cl_error |= clSetKernelArg(pruneKernel, arg_no++, shmem, nullptr);
807     assert(cl_error == CL_SUCCESS);
808
809     cl_error = clEnqueueNDRangeKernel(stream, pruneKernel, 3,
810                                       nullptr, global_work_size, local_work_size,
811                                       0, nullptr, bDoTime ? timer->fetchNextEvent() : nullptr);
812     GMX_RELEASE_ASSERT(CL_SUCCESS == cl_error, ocl_get_error_string(cl_error).c_str());
813
814     if (plist->haveFreshList)
815     {
816         plist->haveFreshList         = false;
817         /* Mark that pruning has been done */
818         nb->timers->didPrune[iloc] = true;
819     }
820     else
821     {
822         /* Mark that rolling pruning has been done */
823         nb->timers->didRollingPrune[iloc] = true;
824     }
825
826     if (bDoTime)
827     {
828         timer->closeTimingRegion(stream);
829     }
830 }
831
832 /*! \brief
833  * Launch asynchronously the download of nonbonded forces from the GPU
834  * (and energies/shift forces if required).
835  */
836 void nbnxn_gpu_launch_cpyback(gmx_nbnxn_ocl_t               *nb,
837                               const struct nbnxn_atomdata_t *nbatom,
838                               int                            flags,
839                               int                            aloc)
840 {
841     cl_int gmx_unused cl_error;
842     int               adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
843
844     /* determine interaction locality from atom locality */
845     int              iloc = gpuAtomToInteractionLocality(aloc);
846
847     cl_atomdata_t   *adat    = nb->atdat;
848     cl_timers_t     *t       = nb->timers;
849     bool             bDoTime = nb->bDoTime;
850     cl_command_queue stream  = nb->stream[iloc];
851
852     bool             bCalcEner   = flags & GMX_FORCE_ENERGY;
853     int              bCalcFshift = flags & GMX_FORCE_VIRIAL;
854
855
856     /* don't launch non-local copy-back if there was no non-local work to do */
857     if (canSkipWork(nb, iloc))
858     {
859         /* TODO An alternative way to signal that non-local work is
860            complete is to use a clEnqueueMarker+clEnqueueBarrier
861            pair. However, the use of bNonLocalStreamActive has the
862            advantage of being local to the host, so probably minimizes
863            overhead. Curiously, for NVIDIA OpenCL with an empty-domain
864            test case, overall simulation performance was higher with
865            the API calls, but this has not been tested on AMD OpenCL,
866            so could be worth considering in future. */
867         nb->bNonLocalStreamActive = false;
868         return;
869     }
870
871     getGpuAtomRange(adat, aloc, adat_begin, adat_len);
872
873     /* beginning of timed D2H section */
874     if (bDoTime)
875     {
876         t->nb_d2h[iloc].openTimingRegion(stream);
877     }
878
879     /* With DD the local D2H transfer can only start after the non-local
880        has been launched. */
881     if (iloc == eintLocal && nb->bNonLocalStreamActive)
882     {
883         sync_ocl_event(stream, &(nb->nonlocal_done));
884     }
885
886     /* DtoH f */
887     ocl_copy_D2H_async(nbatom->out[0].f + adat_begin * 3, adat->f, adat_begin*3*sizeof(float),
888                        (adat_len)* adat->f_elem_size, stream, bDoTime ? t->nb_d2h[iloc].fetchNextEvent() : nullptr);
889
890     /* kick off work */
891     cl_error = clFlush(stream);
892     assert(CL_SUCCESS == cl_error);
893
894     /* After the non-local D2H is launched the nonlocal_done event can be
895        recorded which signals that the local D2H can proceed. This event is not
896        placed after the non-local kernel because we first need the non-local
897        data back first. */
898     if (iloc == eintNonlocal)
899     {
900 #ifdef CL_VERSION_1_2
901         cl_error = clEnqueueMarkerWithWaitList(stream, 0, NULL, &(nb->nonlocal_done));
902 #else
903         cl_error = clEnqueueMarker(stream, &(nb->nonlocal_done));
904 #endif
905         assert(CL_SUCCESS == cl_error);
906         nb->bNonLocalStreamActive = true;
907     }
908
909     /* only transfer energies in the local stream */
910     if (LOCAL_I(iloc))
911     {
912         /* DtoH fshift */
913         if (bCalcFshift)
914         {
915             ocl_copy_D2H_async(nb->nbst.fshift, adat->fshift, 0,
916                                SHIFTS * adat->fshift_elem_size, stream, bDoTime ? t->nb_d2h[iloc].fetchNextEvent() : nullptr);
917         }
918
919         /* DtoH energies */
920         if (bCalcEner)
921         {
922             ocl_copy_D2H_async(nb->nbst.e_lj, adat->e_lj, 0,
923                                sizeof(float), stream, bDoTime ? t->nb_d2h[iloc].fetchNextEvent() : nullptr);
924
925             ocl_copy_D2H_async(nb->nbst.e_el, adat->e_el, 0,
926                                sizeof(float), stream, bDoTime ? t->nb_d2h[iloc].fetchNextEvent() : nullptr);
927         }
928     }
929
930     if (bDoTime)
931     {
932         t->nb_d2h[iloc].closeTimingRegion(stream);
933     }
934 }
935
936 /*! \brief Count pruning kernel time if either kernel has been triggered
937  *
938  *  We do the accounting for either of the two pruning kernel flavors:
939  *   - 1st pass prune: ran during the current step (prior to the force kernel);
940  *   - rolling prune:  ran at the end of the previous step (prior to the current step H2D xq);
941  *
942  * Note that the resetting of cu_timers_t::didPrune and cu_timers_t::didRollingPrune should happen
943  * after calling this function.
944  *
945  * \param[inout] timers   structs with OCL timer objects
946  * \param[inout] timings  GPU task timing data
947  * \param[in] iloc        interaction locality
948  */
949 static void countPruneKernelTime(cl_timers_t               *timers,
950                                  gmx_wallclock_gpu_nbnxn_t *timings,
951                                  const int                  iloc)
952 {
953     // We might have not done any pruning (e.g. if we skipped with empty domains).
954     if (!timers->didPrune[iloc] && !timers->didRollingPrune[iloc])
955     {
956         return;
957     }
958
959     if (timers->didPrune[iloc])
960     {
961         timings->pruneTime.c++;
962         timings->pruneTime.t += timers->prune_k[iloc].getLastRangeTime();
963     }
964
965     if (timers->didRollingPrune[iloc])
966     {
967         timings->dynamicPruneTime.c++;
968         timings->dynamicPruneTime.t += timers->rollingPrune_k[iloc].getLastRangeTime();
969     }
970 }
971
972 /*! \brief Count pruning kernel time if either kernel has been triggered */
973 static void nbnxn_gpu_reduce_staged_outputs(const cl_nb_staging   &nbst,
974                                             int                    iLocality,
975                                             bool                   reduceEnergies,
976                                             bool                   reduceFshift,
977                                             real                  *e_lj,
978                                             real                  *e_el,
979                                             rvec                  *fshift)
980 {
981     /* turn energy calculation always on/off (for debugging/testing only) */
982     reduceEnergies = (reduceEnergies || always_ener) && !never_ener;
983
984     /* add up energies and shift forces (only once at local F wait) */
985     if (LOCAL_I(iLocality))
986     {
987         if (reduceEnergies)
988         {
989             *e_lj += *nbst.e_lj;
990             *e_el += *nbst.e_el;
991         }
992
993         if (reduceFshift)
994         {
995             for (int i = 0; i < SHIFTS; i++)
996             {
997                 fshift[i][0] += (nbst.fshift)[i][0];
998                 fshift[i][1] += (nbst.fshift)[i][1];
999                 fshift[i][2] += (nbst.fshift)[i][2];
1000             }
1001         }
1002     }
1003 }
1004
1005 /*! \brief Do the per-step timing accounting of the nonbonded tasks. */
1006 static void nbnxn_gpu_accumulate_timings(struct gmx_wallclock_gpu_nbnxn_t *timings,
1007                                          cl_timers_t                      *timers,
1008                                          const cl_plist_t                 *plist,
1009                                          int                               atomLocality,
1010                                          bool                              didEnergyKernels,
1011                                          bool                              doTiming)
1012 {
1013     /* timing data accumulation */
1014     if (!doTiming)
1015     {
1016         return;
1017     }
1018
1019     /* determine interaction locality from atom locality */
1020     int iLocality = gpuAtomToInteractionLocality(atomLocality);
1021
1022     /* turn energy calculation always on/off (for debugging/testing only) */
1023     didEnergyKernels = (didEnergyKernels || always_ener) && !never_ener;
1024
1025     /* only increase counter once (at local F wait) */
1026     if (LOCAL_I(iLocality))
1027     {
1028         timings->nb_c++;
1029         timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].c += 1;
1030     }
1031
1032     /* kernel timings */
1033     timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].t +=
1034         timers->nb_k[iLocality].getLastRangeTime();
1035
1036     /* X/q H2D and F D2H timings */
1037     timings->nb_h2d_t += timers->nb_h2d[iLocality].getLastRangeTime();
1038     timings->nb_d2h_t += timers->nb_d2h[iLocality].getLastRangeTime();
1039
1040     /* Count the pruning kernel times for both cases:1st pass (at search step)
1041        and rolling pruning (if called at the previous step).
1042        We do the accounting here as this is the only sync point where we
1043        know (without checking or additional sync-ing) that prune tasks in
1044        in the current stream have completed (having just blocking-waited
1045        for the force D2H). */
1046     countPruneKernelTime(timers, timings, iLocality);
1047
1048     /* only count atdat and pair-list H2D at pair-search step */
1049     if (timers->didPairlistH2D[iLocality])
1050     {
1051         /* atdat transfer timing (add only once, at local F wait) */
1052         if (LOCAL_A(atomLocality))
1053         {
1054             timings->pl_h2d_c++;
1055             timings->pl_h2d_t += timers->atdat.getLastRangeTime();
1056         }
1057
1058         timings->pl_h2d_t += timers->pl_h2d[iLocality].getLastRangeTime();
1059
1060         /* Clear the timing flag for the next step */
1061         timers->didPairlistH2D[iLocality] = false;
1062     }
1063
1064
1065 }
1066
1067 /*! \brief
1068  * Wait for the asynchronously launched nonbonded calculations and data
1069  * transfers to finish.
1070  */
1071 void nbnxn_gpu_wait_for_gpu(gmx_nbnxn_ocl_t *nb,
1072                             int flags, int aloc,
1073                             real *e_lj, real *e_el, rvec *fshift)
1074 {
1075     /* determine interaction locality from atom locality */
1076     int iLocality = gpuAtomToInteractionLocality(aloc);
1077
1078     /* Launch wait/update timers & counters and do reduction into staging buffers
1079        BUT skip it when during the non-local phase there was actually no work to do.
1080        This is consistent with nbnxn_gpu_launch_kernel.
1081
1082        NOTE: if timing with multiple GPUs (streams) becomes possible, the
1083        counters could end up being inconsistent due to not being incremented
1084        on some of the nodes! */
1085     if (!canSkipWork(nb, iLocality))
1086     {
1087         /* Actual sync point. Waits for everything to be finished in the command queue. TODO: Find out if a more fine grained solution is needed */
1088         cl_int gmx_unused cl_error = clFinish(nb->stream[iLocality]);
1089         assert(CL_SUCCESS == cl_error);
1090
1091         bool calcEner   = flags & GMX_FORCE_ENERGY;
1092         bool calcFshift = flags & GMX_FORCE_VIRIAL;
1093
1094         nbnxn_gpu_accumulate_timings(nb->timings, nb->timers, nb->plist[iLocality], aloc, calcEner, nb->bDoTime);
1095
1096         nbnxn_gpu_reduce_staged_outputs(nb->nbst, iLocality, calcEner, calcFshift, e_lj, e_el, fshift);
1097
1098     }
1099
1100     /* Always reset both pruning flags (doesn't hurt doing it even when timing is off). */
1101     nb->timers->didPrune[iLocality] = nb->timers->didRollingPrune[iLocality] = false;
1102
1103     /* Turn off initial list pruning (doesn't hurt if this is not pair-search step). */
1104     nb->plist[iLocality]->haveFreshList = false;
1105 }
1106
1107 /*! \brief Selects the Ewald kernel type, analytical or tabulated, single or twin cut-off. */
1108 int nbnxn_gpu_pick_ewald_kernel_type(bool bTwinCut)
1109 {
1110     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
1111     int  kernel_type;
1112
1113     /* Benchmarking/development environment variables to force the use of
1114        analytical or tabulated Ewald kernel. */
1115     bForceAnalyticalEwald = (getenv("GMX_OCL_NB_ANA_EWALD") != NULL);
1116     bForceTabulatedEwald  = (getenv("GMX_OCL_NB_TAB_EWALD") != NULL);
1117
1118     if (bForceAnalyticalEwald && bForceTabulatedEwald)
1119     {
1120         gmx_incons("Both analytical and tabulated Ewald OpenCL non-bonded kernels "
1121                    "requested through environment variables.");
1122     }
1123
1124     /* OpenCL: By default, use analytical Ewald
1125      * TODO: tabulated does not work, it needs fixing, see init_nbparam() in nbnxn_ocl_data_mgmt.cpp
1126      *
1127      * TODO: decide if dev_info parameter should be added to recognize NVIDIA CC>=3.0 devices.
1128      *
1129      */
1130     //if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
1131     if ((1                         || bForceAnalyticalEwald) && !bForceTabulatedEwald)
1132     {
1133         bUseAnalyticalEwald = true;
1134
1135         if (debug)
1136         {
1137             fprintf(debug, "Using analytical Ewald OpenCL kernels\n");
1138         }
1139     }
1140     else
1141     {
1142         bUseAnalyticalEwald = false;
1143
1144         if (debug)
1145         {
1146             fprintf(debug, "Using tabulated Ewald OpenCL kernels\n");
1147         }
1148     }
1149
1150     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
1151        forces it (use it for debugging/benchmarking only). */
1152     if (!bTwinCut && (getenv("GMX_OCL_NB_EWALD_TWINCUT") == NULL))
1153     {
1154         kernel_type = bUseAnalyticalEwald ? eelOclEWALD_ANA : eelOclEWALD_TAB;
1155     }
1156     else
1157     {
1158         kernel_type = bUseAnalyticalEwald ? eelOclEWALD_ANA_TWIN : eelOclEWALD_TAB_TWIN;
1159     }
1160
1161     return kernel_type;
1162 }