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