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