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