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