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