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