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