cb42ce9ccafdf8bc3ce7a7bd94b0da09ac385311
[alexxy/gromacs.git] / src / gromacs / nbnxm / opencl / nbnxm_ocl_data_mgmt.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_data_mgmt.h
38  *
39  *  \author Anca Hamuraru <anca@streamcomputing.eu>
40  *  \author Dimitrios Karkoulis <dimitris.karkoulis@gmail.com>
41  *  \author Teemu Virolainen <teemu@streamcomputing.eu>
42  *  \author Szilárd Páll <pall.szilard@gmail.com>
43  *  \ingroup module_nbnxm
44  */
45 #include "gmxpre.h"
46
47 #include <assert.h>
48 #include <stdarg.h>
49 #include <stdio.h>
50 #include <stdlib.h>
51 #include <string.h>
52
53 #include <cmath>
54
55 // TODO We would like to move this down, but the way gmx_nbnxn_gpu_t
56 //      is currently declared means this has to be before gpu_types.h
57 #include "nbnxm_ocl_types.h"
58
59 // TODO Remove this comment when the above order issue is resolved
60 #include "gromacs/gpu_utils/gpu_utils.h"
61 #include "gromacs/gpu_utils/oclutils.h"
62 #include "gromacs/hardware/gpu_hw_info.h"
63 #include "gromacs/math/vectypes.h"
64 #include "gromacs/mdlib/force_flags.h"
65 #include "gromacs/mdtypes/interaction_const.h"
66 #include "gromacs/mdtypes/md_enums.h"
67 #include "gromacs/nbnxm/atomdata.h"
68 #include "gromacs/nbnxm/gpu_data_mgmt.h"
69 #include "gromacs/nbnxm/gpu_jit_support.h"
70 #include "gromacs/nbnxm/nbnxm.h"
71 #include "gromacs/nbnxm/nbnxm_gpu.h"
72 #include "gromacs/nbnxm/pairlistsets.h"
73 #include "gromacs/pbcutil/ishift.h"
74 #include "gromacs/timing/gpu_timing.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/real.h"
79 #include "gromacs/utility/smalloc.h"
80
81 #include "nbnxm_ocl_internal.h"
82
83 namespace Nbnxm
84 {
85
86 /*! \brief Copies of values from cl_driver_diagnostics_intel.h,
87  * which isn't guaranteed to be available. */
88 /**@{*/
89 #define CL_CONTEXT_SHOW_DIAGNOSTICS_INTEL 0x4106
90 #define CL_CONTEXT_DIAGNOSTICS_LEVEL_GOOD_INTEL 0x1
91 #define CL_CONTEXT_DIAGNOSTICS_LEVEL_BAD_INTEL 0x2
92 #define CL_CONTEXT_DIAGNOSTICS_LEVEL_NEUTRAL_INTEL 0x4
93 /**@}*/
94
95 /*! \brief This parameter should be determined heuristically from the
96  * kernel execution times
97  *
98  * This value is best for small systems on a single AMD Radeon R9 290X
99  * (and about 5% faster than 40, which is the default for CUDA
100  * devices). Larger simulation systems were quite insensitive to the
101  * value of this parameter.
102  */
103 static unsigned int gpu_min_ci_balanced_factor = 50;
104
105
106 /*! \brief Returns true if LJ combination rules are used in the non-bonded kernels.
107  *
108  * Full doc in nbnxn_ocl_internal.h */
109 bool useLjCombRule(int vdwType)
110 {
111     return (vdwType == evdwOclCUTCOMBGEOM || vdwType == evdwOclCUTCOMBLB);
112 }
113
114 /*! \brief Tabulates the Ewald Coulomb force and initializes the size/scale
115  * and the table GPU array.
116  *
117  * If called with an already allocated table, it just re-uploads the
118  * table.
119  */
120 static void init_ewald_coulomb_force_table(const EwaldCorrectionTables&     tables,
121                                            cl_nbparam_t*                    nbp,
122                                            const gmx_device_runtime_data_t* runData)
123 {
124     cl_mem coul_tab;
125
126     cl_int cl_error;
127
128     if (nbp->coulomb_tab_climg2d != nullptr)
129     {
130         freeDeviceBuffer(&(nbp->coulomb_tab_climg2d));
131     }
132
133     /* Switched from using textures to using buffers */
134     // TODO: decide which alternative is most efficient - textures or buffers.
135     /*
136        cl_image_format array_format;
137
138        array_format.image_channel_data_type = CL_FLOAT;
139        array_format.image_channel_order     = CL_R;
140
141        coul_tab = clCreateImage2D(runData->context, CL_MEM_READ_WRITE | CL_MEM_COPY_HOST_PTR,
142        &array_format, tabsize, 1, 0, ftmp, &cl_error);
143      */
144
145     coul_tab = clCreateBuffer(
146             runData->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY | CL_MEM_COPY_HOST_PTR,
147             tables.tableF.size() * sizeof(cl_float), const_cast<real*>(tables.tableF.data()), &cl_error);
148     GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
149                        ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
150
151     nbp->coulomb_tab_climg2d = coul_tab;
152     nbp->coulomb_tab_scale   = tables.scale;
153 }
154
155
156 /*! \brief Initializes the atomdata structure first time, it only gets filled at
157     pair-search.
158  */
159 static void init_atomdata_first(cl_atomdata_t* ad, int ntypes, gmx_device_runtime_data_t* runData)
160 {
161     cl_int cl_error;
162
163     ad->ntypes = ntypes;
164
165     /* An element of the shift_vec device buffer has the same size as one element
166        of the host side shift_vec buffer. */
167     ad->shift_vec_elem_size = sizeof(*nbnxn_atomdata_t::shift_vec.data());
168
169     ad->shift_vec = clCreateBuffer(runData->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
170                                    SHIFTS * ad->shift_vec_elem_size, nullptr, &cl_error);
171     GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
172                        ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
173     ad->bShiftVecUploaded = CL_FALSE;
174
175     /* An element of the fshift device buffer has the same size as one element
176        of the host side fshift buffer. */
177     ad->fshift_elem_size = sizeof(*cl_nb_staging_t::fshift);
178
179     ad->fshift = clCreateBuffer(runData->context, CL_MEM_READ_WRITE | CL_MEM_HOST_READ_ONLY,
180                                 SHIFTS * ad->fshift_elem_size, nullptr, &cl_error);
181     GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
182                        ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
183
184     ad->e_lj = clCreateBuffer(runData->context, CL_MEM_READ_WRITE | CL_MEM_HOST_READ_ONLY,
185                               sizeof(float), nullptr, &cl_error);
186     GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
187                        ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
188
189     ad->e_el = clCreateBuffer(runData->context, CL_MEM_READ_WRITE | CL_MEM_HOST_READ_ONLY,
190                               sizeof(float), nullptr, &cl_error);
191     GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
192                        ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
193
194     /* initialize to nullptr pointers to data that is not allocated here and will
195        need reallocation in nbnxn_gpu_init_atomdata */
196     ad->xq = nullptr;
197     ad->f  = nullptr;
198
199     /* size -1 indicates that the respective array hasn't been initialized yet */
200     ad->natoms = -1;
201     ad->nalloc = -1;
202 }
203
204 /*! \brief Copies all parameters related to the cut-off from ic to nbp
205  */
206 static void set_cutoff_parameters(cl_nbparam_t* nbp, const interaction_const_t* ic, const PairlistParams& listParams)
207 {
208     nbp->ewald_beta        = ic->ewaldcoeff_q;
209     nbp->sh_ewald          = ic->sh_ewald;
210     nbp->epsfac            = ic->epsfac;
211     nbp->two_k_rf          = 2.0 * ic->k_rf;
212     nbp->c_rf              = ic->c_rf;
213     nbp->rvdw_sq           = ic->rvdw * ic->rvdw;
214     nbp->rcoulomb_sq       = ic->rcoulomb * ic->rcoulomb;
215     nbp->rlistOuter_sq     = listParams.rlistOuter * listParams.rlistOuter;
216     nbp->rlistInner_sq     = listParams.rlistInner * listParams.rlistInner;
217     nbp->useDynamicPruning = listParams.useDynamicPruning;
218
219     nbp->sh_lj_ewald   = ic->sh_lj_ewald;
220     nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
221
222     nbp->rvdw_switch      = ic->rvdw_switch;
223     nbp->dispersion_shift = ic->dispersion_shift;
224     nbp->repulsion_shift  = ic->repulsion_shift;
225     nbp->vdw_switch       = ic->vdw_switch;
226 }
227
228 /*! \brief Returns the kinds of electrostatics and Vdw OpenCL
229  *  kernels that will be used.
230  *
231  * Respectively, these values are from enum eelOcl and enum
232  * evdwOcl. */
233 static void map_interaction_types_to_gpu_kernel_flavors(const interaction_const_t* ic,
234                                                         int                        combRule,
235                                                         int*                       gpu_eeltype,
236                                                         int*                       gpu_vdwtype)
237 {
238     if (ic->vdwtype == evdwCUT)
239     {
240         switch (ic->vdw_modifier)
241         {
242             case eintmodNONE:
243             case eintmodPOTSHIFT:
244                 switch (combRule)
245                 {
246                     case ljcrNONE: *gpu_vdwtype = evdwOclCUT; break;
247                     case ljcrGEOM: *gpu_vdwtype = evdwOclCUTCOMBGEOM; break;
248                     case ljcrLB: *gpu_vdwtype = evdwOclCUTCOMBLB; break;
249                     default:
250                         gmx_incons(
251                                 "The requested LJ combination rule is not implemented in the "
252                                 "OpenCL GPU accelerated kernels!");
253                 }
254                 break;
255             case eintmodFORCESWITCH: *gpu_vdwtype = evdwOclFSWITCH; break;
256             case eintmodPOTSWITCH: *gpu_vdwtype = evdwOclPSWITCH; break;
257             default:
258                 gmx_incons(
259                         "The requested VdW interaction modifier is not implemented in the GPU "
260                         "accelerated kernels!");
261         }
262     }
263     else if (ic->vdwtype == evdwPME)
264     {
265         if (ic->ljpme_comb_rule == ljcrGEOM)
266         {
267             *gpu_vdwtype = evdwOclEWALDGEOM;
268         }
269         else
270         {
271             *gpu_vdwtype = evdwOclEWALDLB;
272         }
273     }
274     else
275     {
276         gmx_incons("The requested VdW type is not implemented in the GPU accelerated kernels!");
277     }
278
279     if (ic->eeltype == eelCUT)
280     {
281         *gpu_eeltype = eelOclCUT;
282     }
283     else if (EEL_RF(ic->eeltype))
284     {
285         *gpu_eeltype = eelOclRF;
286     }
287     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
288     {
289         *gpu_eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
290     }
291     else
292     {
293         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
294         gmx_incons(
295                 "The requested electrostatics type is not implemented in the GPU accelerated "
296                 "kernels!");
297     }
298 }
299
300 /*! \brief Initializes the nonbonded parameter data structure.
301  */
302 static void init_nbparam(cl_nbparam_t*                    nbp,
303                          const interaction_const_t*       ic,
304                          const PairlistParams&            listParams,
305                          const nbnxn_atomdata_t::Params&  nbatParams,
306                          const gmx_device_runtime_data_t* runData)
307 {
308     cl_int cl_error;
309
310     set_cutoff_parameters(nbp, ic, listParams);
311
312     map_interaction_types_to_gpu_kernel_flavors(ic, nbatParams.comb_rule, &(nbp->eeltype), &(nbp->vdwtype));
313
314     if (ic->vdwtype == evdwPME)
315     {
316         if (ic->ljpme_comb_rule == ljcrGEOM)
317         {
318             GMX_ASSERT(nbatParams.comb_rule == ljcrGEOM, "Combination rule mismatch!");
319         }
320         else
321         {
322             GMX_ASSERT(nbatParams.comb_rule == ljcrLB, "Combination rule mismatch!");
323         }
324     }
325     /* generate table for PME */
326     nbp->coulomb_tab_climg2d = nullptr;
327     if (nbp->eeltype == eelOclEWALD_TAB || nbp->eeltype == eelOclEWALD_TAB_TWIN)
328     {
329         GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
330         init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, runData);
331     }
332     else
333     // TODO: improvement needed.
334     // The image2d is created here even if eeltype is not eelCuEWALD_TAB or eelCuEWALD_TAB_TWIN
335     // because the OpenCL kernels don't accept nullptr values for image2D parameters.
336     {
337         /* Switched from using textures to using buffers */
338         // TODO: decide which alternative is most efficient - textures or buffers.
339         /*
340            cl_image_format array_format;
341
342            array_format.image_channel_data_type = CL_FLOAT;
343            array_format.image_channel_order     = CL_R;
344
345            nbp->coulomb_tab_climg2d = clCreateImage2D(runData->context, CL_MEM_READ_WRITE,
346             &array_format, 1, 1, 0, nullptr, &cl_error);
347          */
348
349         nbp->coulomb_tab_climg2d = clCreateBuffer(runData->context, CL_MEM_READ_ONLY,
350                                                   sizeof(cl_float), nullptr, &cl_error);
351         GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
352                            ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
353     }
354
355     const int nnbfp      = 2 * nbatParams.numTypes * nbatParams.numTypes;
356     const int nnbfp_comb = 2 * nbatParams.numTypes;
357
358     {
359         /* Switched from using textures to using buffers */
360         // TODO: decide which alternative is most efficient - textures or buffers.
361         /*
362            cl_image_format array_format;
363
364            array_format.image_channel_data_type = CL_FLOAT;
365            array_format.image_channel_order     = CL_R;
366
367            nbp->nbfp_climg2d = clCreateImage2D(runData->context, CL_MEM_READ_ONLY |
368            CL_MEM_COPY_HOST_PTR, &array_format, nnbfp, 1, 0, nbat->nbfp, &cl_error);
369          */
370
371         nbp->nbfp_climg2d = clCreateBuffer(
372                 runData->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY | CL_MEM_COPY_HOST_PTR,
373                 nnbfp * sizeof(cl_float), const_cast<float*>(nbatParams.nbfp.data()), &cl_error);
374         GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
375                            ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
376
377         if (ic->vdwtype == evdwPME)
378         {
379             /* Switched from using textures to using buffers */
380             // TODO: decide which alternative is most efficient - textures or buffers.
381             /*  nbp->nbfp_comb_climg2d = clCreateImage2D(runData->context, CL_MEM_READ_WRITE |
382                CL_MEM_COPY_HOST_PTR, &array_format, nnbfp_comb, 1, 0, nbat->nbfp_comb, &cl_error);*/
383             nbp->nbfp_comb_climg2d = clCreateBuffer(
384                     runData->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY | CL_MEM_COPY_HOST_PTR,
385                     nnbfp_comb * sizeof(cl_float), const_cast<float*>(nbatParams.nbfp_comb.data()),
386                     &cl_error);
387             GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
388                                ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
389         }
390         else
391         {
392             // TODO: improvement needed.
393             // The image2d is created here even if vdwtype is not evdwPME because the OpenCL kernels
394             // don't accept nullptr values for image2D parameters.
395             /* Switched from using textures to using buffers */
396             // TODO: decide which alternative is most efficient - textures or buffers.
397             /* nbp->nbfp_comb_climg2d = clCreateImage2D(runData->context, CL_MEM_READ_WRITE,
398                 &array_format, 1, 1, 0, nullptr, &cl_error);*/
399             nbp->nbfp_comb_climg2d = clCreateBuffer(runData->context, CL_MEM_READ_ONLY,
400                                                     sizeof(cl_float), nullptr, &cl_error);
401             GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
402                                ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
403         }
404     }
405 }
406
407 //! This function is documented in the header file
408 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t* nbv, const interaction_const_t* ic)
409 {
410     if (!nbv || !nbv->useGpu())
411     {
412         return;
413     }
414     gmx_nbnxn_ocl_t* nb  = nbv->gpu_nbv;
415     cl_nbparam_t*    nbp = nb->nbparam;
416
417     set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
418
419     nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
420
421     GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
422     init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, nb->dev_rundata);
423 }
424
425 /*! \brief Initializes the pair list data structure.
426  */
427 static void init_plist(cl_plist_t* pl)
428 {
429     /* initialize to nullptr pointers to data that is not allocated here and will
430        need reallocation in nbnxn_gpu_init_pairlist */
431     pl->sci   = nullptr;
432     pl->cj4   = nullptr;
433     pl->imask = nullptr;
434     pl->excl  = nullptr;
435
436     /* size -1 indicates that the respective array hasn't been initialized yet */
437     pl->na_c          = -1;
438     pl->nsci          = -1;
439     pl->sci_nalloc    = -1;
440     pl->ncj4          = -1;
441     pl->cj4_nalloc    = -1;
442     pl->nimask        = -1;
443     pl->imask_nalloc  = -1;
444     pl->nexcl         = -1;
445     pl->excl_nalloc   = -1;
446     pl->haveFreshList = false;
447 }
448
449 /*! \brief Initializes the timings data structure.
450  */
451 static void init_timings(gmx_wallclock_gpu_nbnxn_t* t)
452 {
453     int i, j;
454
455     t->nb_h2d_t = 0.0;
456     t->nb_d2h_t = 0.0;
457     t->nb_c     = 0;
458     t->pl_h2d_t = 0.0;
459     t->pl_h2d_c = 0;
460     for (i = 0; i < 2; i++)
461     {
462         for (j = 0; j < 2; j++)
463         {
464             t->ktime[i][j].t = 0.0;
465             t->ktime[i][j].c = 0;
466         }
467     }
468
469     t->pruneTime.c        = 0;
470     t->pruneTime.t        = 0.0;
471     t->dynamicPruneTime.c = 0;
472     t->dynamicPruneTime.t = 0.0;
473 }
474
475
476 //! OpenCL notification callback function
477 static void CL_CALLBACK ocl_notify_fn(const char* pErrInfo,
478                                       const void gmx_unused* private_info,
479                                       size_t gmx_unused cb,
480                                       void gmx_unused* user_data)
481 {
482     if (pErrInfo != nullptr)
483     {
484         printf("%s\n", pErrInfo); // Print error/hint
485     }
486 }
487
488 /*! \brief Creates context for OpenCL GPU given by \p mygpu
489  *
490  * A fatal error results if creation fails.
491  *
492  * \param[inout] runtimeData runtime data including program and context
493  * \param[in]    devInfo     device info struct
494  * \param[in]    rank        MPI rank (for error reporting)
495  */
496 static void nbnxn_gpu_create_context(gmx_device_runtime_data_t* runtimeData,
497                                      const gmx_device_info_t*   devInfo,
498                                      int                        rank)
499 {
500     cl_context_properties context_properties[5];
501     cl_platform_id        platform_id;
502     cl_device_id          device_id;
503     cl_context            context;
504     cl_int                cl_error;
505
506     assert(runtimeData != nullptr);
507     assert(devInfo != nullptr);
508
509     platform_id = devInfo->ocl_gpu_id.ocl_platform_id;
510     device_id   = devInfo->ocl_gpu_id.ocl_device_id;
511
512     int i                   = 0;
513     context_properties[i++] = CL_CONTEXT_PLATFORM;
514     context_properties[i++] = reinterpret_cast<cl_context_properties>(platform_id);
515     if (getenv("GMX_OCL_SHOW_DIAGNOSTICS"))
516     {
517         context_properties[i++] = CL_CONTEXT_SHOW_DIAGNOSTICS_INTEL;
518         context_properties[i++] =
519                 CL_CONTEXT_DIAGNOSTICS_LEVEL_BAD_INTEL | CL_CONTEXT_DIAGNOSTICS_LEVEL_NEUTRAL_INTEL;
520     }
521     context_properties[i++] = 0; /* Terminates the list of properties */
522
523     context = clCreateContext(context_properties, 1, &device_id, ocl_notify_fn, nullptr, &cl_error);
524     if (CL_SUCCESS != cl_error)
525     {
526         gmx_fatal(FARGS, "On rank %d failed to create context for GPU #%s:\n OpenCL error %d: %s",
527                   rank, devInfo->device_name, cl_error, ocl_get_error_string(cl_error).c_str());
528     }
529
530     runtimeData->context = context;
531 }
532
533 /*! \brief Initializes the OpenCL kernel pointers of the nbnxn_ocl_ptr_t input data structure. */
534 static cl_kernel nbnxn_gpu_create_kernel(gmx_nbnxn_ocl_t* nb, const char* kernel_name)
535 {
536     cl_kernel kernel;
537     cl_int    cl_error;
538
539     kernel = clCreateKernel(nb->dev_rundata->program, kernel_name, &cl_error);
540     if (CL_SUCCESS != cl_error)
541     {
542         gmx_fatal(FARGS, "Failed to create kernel '%s' for GPU #%s: OpenCL error %d", kernel_name,
543                   nb->dev_info->device_name, cl_error);
544     }
545
546     return kernel;
547 }
548
549 /*! \brief Clears nonbonded shift force output array and energy outputs on the GPU.
550  */
551 static void nbnxn_ocl_clear_e_fshift(gmx_nbnxn_ocl_t* nb)
552 {
553
554     cl_int           cl_error;
555     cl_atomdata_t*   adat = nb->atdat;
556     cl_command_queue ls   = nb->stream[InteractionLocality::Local];
557
558     size_t local_work_size[3]  = { 1, 1, 1 };
559     size_t global_work_size[3] = { 1, 1, 1 };
560
561     cl_int shifts = SHIFTS * 3;
562
563     cl_int arg_no;
564
565     cl_kernel zero_e_fshift = nb->kernel_zero_e_fshift;
566
567     local_work_size[0] = 64;
568     // Round the total number of threads up from the array size
569     global_work_size[0] = ((shifts + local_work_size[0] - 1) / local_work_size[0]) * local_work_size[0];
570
571     arg_no   = 0;
572     cl_error = clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_mem), &(adat->fshift));
573     cl_error |= clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_mem), &(adat->e_lj));
574     cl_error |= clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_mem), &(adat->e_el));
575     cl_error |= clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_uint), &shifts);
576     GMX_ASSERT(cl_error == CL_SUCCESS, ocl_get_error_string(cl_error).c_str());
577
578     cl_error = clEnqueueNDRangeKernel(ls, zero_e_fshift, 3, nullptr, global_work_size,
579                                       local_work_size, 0, nullptr, nullptr);
580     GMX_ASSERT(cl_error == CL_SUCCESS, ocl_get_error_string(cl_error).c_str());
581 }
582
583 /*! \brief Initializes the OpenCL kernel pointers of the nbnxn_ocl_ptr_t input data structure. */
584 static void nbnxn_gpu_init_kernels(gmx_nbnxn_ocl_t* nb)
585 {
586     /* Init to 0 main kernel arrays */
587     /* They will be later on initialized in select_nbnxn_kernel */
588     // TODO: consider always creating all variants of the kernels here so that there is no
589     // need for late call to clCreateKernel -- if that gives any advantage?
590     memset(nb->kernel_ener_noprune_ptr, 0, sizeof(nb->kernel_ener_noprune_ptr));
591     memset(nb->kernel_ener_prune_ptr, 0, sizeof(nb->kernel_ener_prune_ptr));
592     memset(nb->kernel_noener_noprune_ptr, 0, sizeof(nb->kernel_noener_noprune_ptr));
593     memset(nb->kernel_noener_prune_ptr, 0, sizeof(nb->kernel_noener_prune_ptr));
594
595     /* Init pruning kernels
596      *
597      * TODO: we could avoid creating kernels if dynamic pruning is turned off,
598      * but ATM that depends on force flags not passed into the initialization.
599      */
600     nb->kernel_pruneonly[epruneFirst] = nbnxn_gpu_create_kernel(nb, "nbnxn_kernel_prune_opencl");
601     nb->kernel_pruneonly[epruneRolling] =
602             nbnxn_gpu_create_kernel(nb, "nbnxn_kernel_prune_rolling_opencl");
603
604     /* Init auxiliary kernels */
605     nb->kernel_zero_e_fshift = nbnxn_gpu_create_kernel(nb, "zero_e_fshift");
606 }
607
608 /*! \brief Initializes simulation constant data.
609  *
610  *  Initializes members of the atomdata and nbparam structs and
611  *  clears e/fshift output buffers.
612  */
613 static void nbnxn_ocl_init_const(gmx_nbnxn_ocl_t*                nb,
614                                  const interaction_const_t*      ic,
615                                  const PairlistParams&           listParams,
616                                  const nbnxn_atomdata_t::Params& nbatParams)
617 {
618     init_atomdata_first(nb->atdat, nbatParams.numTypes, nb->dev_rundata);
619     init_nbparam(nb->nbparam, ic, listParams, nbatParams, nb->dev_rundata);
620 }
621
622
623 //! This function is documented in the header file
624 gmx_nbnxn_ocl_t* gpu_init(const gmx_device_info_t*   deviceInfo,
625                           const interaction_const_t* ic,
626                           const PairlistParams&      listParams,
627                           const nbnxn_atomdata_t*    nbat,
628                           const int                  rank,
629                           const gmx_bool             bLocalAndNonlocal)
630 {
631     gmx_nbnxn_ocl_t*            nb;
632     cl_int                      cl_error;
633     cl_command_queue_properties queue_properties;
634
635     assert(ic);
636
637     snew(nb, 1);
638     snew(nb->atdat, 1);
639     snew(nb->nbparam, 1);
640     snew(nb->plist[InteractionLocality::Local], 1);
641     if (bLocalAndNonlocal)
642     {
643         snew(nb->plist[InteractionLocality::NonLocal], 1);
644     }
645
646     nb->bUseTwoStreams = static_cast<cl_bool>(bLocalAndNonlocal);
647
648     nb->timers = new cl_timers_t();
649     snew(nb->timings, 1);
650
651     /* set device info, just point it to the right GPU among the detected ones */
652     nb->dev_info = deviceInfo;
653     snew(nb->dev_rundata, 1);
654
655     /* init nbst */
656     pmalloc(reinterpret_cast<void**>(&nb->nbst.e_lj), sizeof(*nb->nbst.e_lj));
657     pmalloc(reinterpret_cast<void**>(&nb->nbst.e_el), sizeof(*nb->nbst.e_el));
658     pmalloc(reinterpret_cast<void**>(&nb->nbst.fshift), SHIFTS * sizeof(*nb->nbst.fshift));
659
660     init_plist(nb->plist[InteractionLocality::Local]);
661
662     /* OpenCL timing disabled if GMX_DISABLE_GPU_TIMING is defined. */
663     nb->bDoTime = static_cast<cl_bool>(getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
664
665     /* Create queues only after bDoTime has been initialized */
666     if (nb->bDoTime)
667     {
668         queue_properties = CL_QUEUE_PROFILING_ENABLE;
669     }
670     else
671     {
672         queue_properties = 0;
673     }
674
675     nbnxn_gpu_create_context(nb->dev_rundata, nb->dev_info, rank);
676
677     /* local/non-local GPU streams */
678     nb->stream[InteractionLocality::Local] = clCreateCommandQueue(
679             nb->dev_rundata->context, nb->dev_info->ocl_gpu_id.ocl_device_id, queue_properties, &cl_error);
680     if (CL_SUCCESS != cl_error)
681     {
682         gmx_fatal(FARGS, "On rank %d failed to create context for GPU #%s: OpenCL error %d", rank,
683                   nb->dev_info->device_name, cl_error);
684     }
685
686     if (nb->bUseTwoStreams)
687     {
688         init_plist(nb->plist[InteractionLocality::NonLocal]);
689
690         nb->stream[InteractionLocality::NonLocal] =
691                 clCreateCommandQueue(nb->dev_rundata->context, nb->dev_info->ocl_gpu_id.ocl_device_id,
692                                      queue_properties, &cl_error);
693         if (CL_SUCCESS != cl_error)
694         {
695             gmx_fatal(FARGS, "On rank %d failed to create context for GPU #%s: OpenCL error %d",
696                       rank, nb->dev_info->device_name, cl_error);
697         }
698     }
699
700     if (nb->bDoTime)
701     {
702         init_timings(nb->timings);
703     }
704
705     nbnxn_ocl_init_const(nb, ic, listParams, nbat->params());
706
707     /* Enable LJ param manual prefetch for AMD or Intel or if we request through env. var.
708      * TODO: decide about NVIDIA
709      */
710     nb->bPrefetchLjParam = (getenv("GMX_OCL_DISABLE_I_PREFETCH") == nullptr)
711                            && ((nb->dev_info->vendor_e == OCL_VENDOR_AMD)
712                                || (nb->dev_info->vendor_e == OCL_VENDOR_INTEL)
713                                || (getenv("GMX_OCL_ENABLE_I_PREFETCH") != nullptr));
714
715     /* NOTE: in CUDA we pick L1 cache configuration for the nbnxn kernels here,
716      * but sadly this is not supported in OpenCL (yet?). Consider adding it if
717      * it becomes supported.
718      */
719     nbnxn_gpu_compile_kernels(nb);
720     nbnxn_gpu_init_kernels(nb);
721
722     /* clear energy and shift force outputs */
723     nbnxn_ocl_clear_e_fshift(nb);
724
725     if (debug)
726     {
727         fprintf(debug, "Initialized OpenCL data structures.\n");
728     }
729
730     return nb;
731 }
732
733 /*! \brief Clears the first natoms_clear elements of the GPU nonbonded force output array.
734  */
735 static void nbnxn_ocl_clear_f(gmx_nbnxn_ocl_t* nb, int natoms_clear)
736 {
737     if (natoms_clear == 0)
738     {
739         return;
740     }
741
742     cl_int gmx_used_in_debug cl_error;
743
744     cl_atomdata_t*   atomData = nb->atdat;
745     cl_command_queue ls       = nb->stream[InteractionLocality::Local];
746     cl_float         value    = 0.0F;
747
748     cl_error = clEnqueueFillBuffer(ls, atomData->f, &value, sizeof(cl_float), 0,
749                                    natoms_clear * sizeof(rvec), 0, nullptr, nullptr);
750     GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
751                        ("clEnqueueFillBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
752 }
753
754 //! This function is documented in the header file
755 void gpu_clear_outputs(gmx_nbnxn_ocl_t* nb, bool computeVirial)
756 {
757     nbnxn_ocl_clear_f(nb, nb->atdat->natoms);
758     /* clear shift force array and energies if the outputs were
759        used in the current step */
760     if (computeVirial)
761     {
762         nbnxn_ocl_clear_e_fshift(nb);
763     }
764
765     /* kick off buffer clearing kernel to ensure concurrency with constraints/update */
766     cl_int gmx_unused cl_error;
767     cl_error = clFlush(nb->stream[InteractionLocality::Local]);
768     assert(CL_SUCCESS == cl_error);
769 }
770
771 //! This function is documented in the header file
772 void gpu_init_pairlist(gmx_nbnxn_ocl_t* nb, const NbnxnPairlistGpu* h_plist, const InteractionLocality iloc)
773 {
774     char sbuf[STRLEN];
775     // Timing accumulation should happen only if there was work to do
776     // because getLastRangeTime() gets skipped with empty lists later
777     // which leads to the counter not being reset.
778     bool             bDoTime = ((nb->bDoTime == CL_TRUE) && !h_plist->sci.empty());
779     cl_command_queue stream  = nb->stream[iloc];
780     cl_plist_t*      d_plist = nb->plist[iloc];
781
782     if (d_plist->na_c < 0)
783     {
784         d_plist->na_c = h_plist->na_ci;
785     }
786     else
787     {
788         if (d_plist->na_c != h_plist->na_ci)
789         {
790             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
791                     d_plist->na_c, h_plist->na_ci);
792             gmx_incons(sbuf);
793         }
794     }
795
796     gpu_timers_t::Interaction& iTimers = nb->timers->interaction[iloc];
797
798     if (bDoTime)
799     {
800         iTimers.pl_h2d.openTimingRegion(stream);
801         iTimers.didPairlistH2D = true;
802     }
803
804     // TODO most of this function is same in CUDA and OpenCL, move into the header
805     DeviceContext context = nb->dev_rundata->context;
806
807     reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(), &d_plist->nsci, &d_plist->sci_nalloc, context);
808     copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(), stream,
809                        GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
810
811     reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(), &d_plist->ncj4, &d_plist->cj4_nalloc, context);
812     copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(), stream,
813                        GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
814
815     reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size() * c_nbnxnGpuClusterpairSplit,
816                            &d_plist->nimask, &d_plist->imask_nalloc, context);
817
818     reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(), &d_plist->nexcl,
819                            &d_plist->excl_nalloc, context);
820     copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(), stream,
821                        GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
822
823     if (bDoTime)
824     {
825         iTimers.pl_h2d.closeTimingRegion(stream);
826     }
827
828     /* need to prune the pair list during the next step */
829     d_plist->haveFreshList = true;
830 }
831
832 //! This function is documented in the header file
833 void gpu_upload_shiftvec(gmx_nbnxn_ocl_t* nb, const nbnxn_atomdata_t* nbatom)
834 {
835     cl_atomdata_t*   adat = nb->atdat;
836     cl_command_queue ls   = nb->stream[InteractionLocality::Local];
837
838     /* only if we have a dynamic box */
839     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
840     {
841         ocl_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(), 0,
842                            SHIFTS * adat->shift_vec_elem_size, ls, nullptr);
843         adat->bShiftVecUploaded = CL_TRUE;
844     }
845 }
846
847 //! This function is documented in the header file
848 void gpu_init_atomdata(gmx_nbnxn_ocl_t* nb, const nbnxn_atomdata_t* nbat)
849 {
850     cl_int           cl_error;
851     int              nalloc, natoms;
852     bool             realloced;
853     bool             bDoTime = nb->bDoTime == CL_TRUE;
854     cl_timers_t*     timers  = nb->timers;
855     cl_atomdata_t*   d_atdat = nb->atdat;
856     cl_command_queue ls      = nb->stream[InteractionLocality::Local];
857
858     natoms    = nbat->numAtoms();
859     realloced = false;
860
861     if (bDoTime)
862     {
863         /* time async copy */
864         timers->atdat.openTimingRegion(ls);
865     }
866
867     /* need to reallocate if we have to copy more atoms than the amount of space
868        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
869     if (natoms > d_atdat->nalloc)
870     {
871         nalloc = over_alloc_small(natoms);
872
873         /* free up first if the arrays have already been initialized */
874         if (d_atdat->nalloc != -1)
875         {
876             freeDeviceBuffer(&d_atdat->f);
877             freeDeviceBuffer(&d_atdat->xq);
878             freeDeviceBuffer(&d_atdat->lj_comb);
879             freeDeviceBuffer(&d_atdat->atom_types);
880         }
881
882         d_atdat->f_elem_size = sizeof(rvec);
883
884         d_atdat->f = clCreateBuffer(nb->dev_rundata->context, CL_MEM_READ_WRITE | CL_MEM_HOST_READ_ONLY,
885                                     nalloc * d_atdat->f_elem_size, nullptr, &cl_error);
886         GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
887                            ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
888
889         d_atdat->xq = clCreateBuffer(nb->dev_rundata->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
890                                      nalloc * sizeof(cl_float4), nullptr, &cl_error);
891         GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
892                            ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
893
894         if (useLjCombRule(nb->nbparam->vdwtype))
895         {
896             d_atdat->lj_comb = clCreateBuffer(nb->dev_rundata->context,
897                                               CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
898                                               nalloc * sizeof(cl_float2), nullptr, &cl_error);
899             GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
900                                ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
901         }
902         else
903         {
904             d_atdat->atom_types = clCreateBuffer(nb->dev_rundata->context,
905                                                  CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
906                                                  nalloc * sizeof(int), nullptr, &cl_error);
907             GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
908                                ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
909         }
910
911         d_atdat->nalloc = nalloc;
912         realloced       = true;
913     }
914
915     d_atdat->natoms       = natoms;
916     d_atdat->natoms_local = nbat->natoms_local;
917
918     /* need to clear GPU f output if realloc happened */
919     if (realloced)
920     {
921         nbnxn_ocl_clear_f(nb, nalloc);
922     }
923
924     if (useLjCombRule(nb->nbparam->vdwtype))
925     {
926         ocl_copy_H2D_async(d_atdat->lj_comb, nbat->params().lj_comb.data(), 0, natoms * sizeof(cl_float2),
927                            ls, bDoTime ? timers->atdat.fetchNextEvent() : nullptr);
928     }
929     else
930     {
931         ocl_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(), 0, natoms * sizeof(int),
932                            ls, bDoTime ? timers->atdat.fetchNextEvent() : nullptr);
933     }
934
935     if (bDoTime)
936     {
937         timers->atdat.closeTimingRegion(ls);
938     }
939
940     /* kick off the tasks enqueued above to ensure concurrency with the search */
941     cl_error = clFlush(ls);
942     GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
943                        ("clFlush failed: " + ocl_get_error_string(cl_error)).c_str());
944 }
945
946 /*! \brief Releases an OpenCL kernel pointer */
947 static void free_kernel(cl_kernel* kernel_ptr)
948 {
949     cl_int gmx_unused cl_error;
950
951     assert(nullptr != kernel_ptr);
952
953     if (*kernel_ptr)
954     {
955         cl_error = clReleaseKernel(*kernel_ptr);
956         GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
957                            ("clReleaseKernel failed: " + ocl_get_error_string(cl_error)).c_str());
958
959         *kernel_ptr = nullptr;
960     }
961 }
962
963 /*! \brief Releases a list of OpenCL kernel pointers */
964 static void free_kernels(cl_kernel* kernels, int count)
965 {
966     int i;
967
968     for (i = 0; i < count; i++)
969     {
970         free_kernel(kernels + i);
971     }
972 }
973
974 /*! \brief Free the OpenCL runtime data (context and program).
975  *
976  *  The function releases the OpenCL context and program assuciated with the
977  *  device that the calling PP rank is running on.
978  *
979  *  \param runData [in]  porinter to the structure with runtime data.
980  */
981 static void free_gpu_device_runtime_data(gmx_device_runtime_data_t* runData)
982 {
983     if (runData == nullptr)
984     {
985         return;
986     }
987
988     cl_int gmx_unused cl_error;
989
990     if (runData->context)
991     {
992         cl_error = clReleaseContext(runData->context);
993         GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
994                            ("clReleaseContext failed: " + ocl_get_error_string(cl_error)).c_str());
995         runData->context = nullptr;
996     }
997
998     if (runData->program)
999     {
1000         cl_error = clReleaseProgram(runData->program);
1001         GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
1002                            ("clReleaseProgram failed: " + ocl_get_error_string(cl_error)).c_str());
1003         runData->program = nullptr;
1004     }
1005 }
1006
1007 //! This function is documented in the header file
1008 void gpu_free(gmx_nbnxn_ocl_t* nb)
1009 {
1010     if (nb == nullptr)
1011     {
1012         return;
1013     }
1014
1015     /* Free kernels */
1016     int kernel_count = sizeof(nb->kernel_ener_noprune_ptr) / sizeof(nb->kernel_ener_noprune_ptr[0][0]);
1017     free_kernels(nb->kernel_ener_noprune_ptr[0], kernel_count);
1018
1019     kernel_count = sizeof(nb->kernel_ener_prune_ptr) / sizeof(nb->kernel_ener_prune_ptr[0][0]);
1020     free_kernels(nb->kernel_ener_prune_ptr[0], kernel_count);
1021
1022     kernel_count = sizeof(nb->kernel_noener_noprune_ptr) / sizeof(nb->kernel_noener_noprune_ptr[0][0]);
1023     free_kernels(nb->kernel_noener_noprune_ptr[0], kernel_count);
1024
1025     kernel_count = sizeof(nb->kernel_noener_prune_ptr) / sizeof(nb->kernel_noener_prune_ptr[0][0]);
1026     free_kernels(nb->kernel_noener_prune_ptr[0], kernel_count);
1027
1028     free_kernel(&(nb->kernel_zero_e_fshift));
1029
1030     /* Free atdat */
1031     freeDeviceBuffer(&(nb->atdat->xq));
1032     freeDeviceBuffer(&(nb->atdat->f));
1033     freeDeviceBuffer(&(nb->atdat->e_lj));
1034     freeDeviceBuffer(&(nb->atdat->e_el));
1035     freeDeviceBuffer(&(nb->atdat->fshift));
1036     freeDeviceBuffer(&(nb->atdat->lj_comb));
1037     freeDeviceBuffer(&(nb->atdat->atom_types));
1038     freeDeviceBuffer(&(nb->atdat->shift_vec));
1039     sfree(nb->atdat);
1040
1041     /* Free nbparam */
1042     freeDeviceBuffer(&(nb->nbparam->nbfp_climg2d));
1043     freeDeviceBuffer(&(nb->nbparam->nbfp_comb_climg2d));
1044     freeDeviceBuffer(&(nb->nbparam->coulomb_tab_climg2d));
1045     sfree(nb->nbparam);
1046
1047     /* Free plist */
1048     auto* plist = nb->plist[InteractionLocality::Local];
1049     freeDeviceBuffer(&plist->sci);
1050     freeDeviceBuffer(&plist->cj4);
1051     freeDeviceBuffer(&plist->imask);
1052     freeDeviceBuffer(&plist->excl);
1053     sfree(plist);
1054     if (nb->bUseTwoStreams)
1055     {
1056         auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
1057         freeDeviceBuffer(&plist_nl->sci);
1058         freeDeviceBuffer(&plist_nl->cj4);
1059         freeDeviceBuffer(&plist_nl->imask);
1060         freeDeviceBuffer(&plist_nl->excl);
1061         sfree(plist_nl);
1062     }
1063
1064     /* Free nbst */
1065     pfree(nb->nbst.e_lj);
1066     nb->nbst.e_lj = nullptr;
1067
1068     pfree(nb->nbst.e_el);
1069     nb->nbst.e_el = nullptr;
1070
1071     pfree(nb->nbst.fshift);
1072     nb->nbst.fshift = nullptr;
1073
1074     /* Free command queues */
1075     clReleaseCommandQueue(nb->stream[InteractionLocality::Local]);
1076     nb->stream[InteractionLocality::Local] = nullptr;
1077     if (nb->bUseTwoStreams)
1078     {
1079         clReleaseCommandQueue(nb->stream[InteractionLocality::NonLocal]);
1080         nb->stream[InteractionLocality::NonLocal] = nullptr;
1081     }
1082     /* Free other events */
1083     if (nb->nonlocal_done)
1084     {
1085         clReleaseEvent(nb->nonlocal_done);
1086         nb->nonlocal_done = nullptr;
1087     }
1088     if (nb->misc_ops_and_local_H2D_done)
1089     {
1090         clReleaseEvent(nb->misc_ops_and_local_H2D_done);
1091         nb->misc_ops_and_local_H2D_done = nullptr;
1092     }
1093
1094     free_gpu_device_runtime_data(nb->dev_rundata);
1095     sfree(nb->dev_rundata);
1096
1097     /* Free timers and timings */
1098     delete nb->timers;
1099     sfree(nb->timings);
1100     sfree(nb);
1101
1102     if (debug)
1103     {
1104         fprintf(debug, "Cleaned up OpenCL data structures.\n");
1105     }
1106 }
1107
1108 //! This function is documented in the header file
1109 gmx_wallclock_gpu_nbnxn_t* gpu_get_timings(gmx_nbnxn_ocl_t* nb)
1110 {
1111     return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
1112 }
1113
1114 //! This function is documented in the header file
1115 void gpu_reset_timings(nonbonded_verlet_t* nbv)
1116 {
1117     if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
1118     {
1119         init_timings(nbv->gpu_nbv->timings);
1120     }
1121 }
1122
1123 //! This function is documented in the header file
1124 int gpu_min_ci_balanced(gmx_nbnxn_ocl_t* nb)
1125 {
1126     return nb != nullptr ? gpu_min_ci_balanced_factor * nb->dev_info->compute_units : 0;
1127 }
1128
1129 //! This function is documented in the header file
1130 gmx_bool gpu_is_kernel_ewald_analytical(const gmx_nbnxn_ocl_t* nb)
1131 {
1132     return ((nb->nbparam->eeltype == eelOclEWALD_ANA) || (nb->nbparam->eeltype == eelOclEWALD_ANA_TWIN));
1133 }
1134
1135 } // namespace Nbnxm