Split lines with many copyright years
[gromacs.git] / src / gromacs / nbnxm / opencl / nbnxm_ocl_data_mgmt.cpp
blobcb42ce9ccafdf8bc3ce7a7bd94b0da09ac385311
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
36 /*! \internal \file
37 * \brief Define OpenCL implementation of nbnxm_gpu_data_mgmt.h
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
45 #include "gmxpre.h"
47 #include <assert.h>
48 #include <stdarg.h>
49 #include <stdio.h>
50 #include <stdlib.h>
51 #include <string.h>
53 #include <cmath>
55 // TODO We would like to move this down, but the way gmx_nbnxn_gpu_t
56 // is currently declared means this has to be before gpu_types.h
57 #include "nbnxm_ocl_types.h"
59 // TODO Remove this comment when the above order issue is resolved
60 #include "gromacs/gpu_utils/gpu_utils.h"
61 #include "gromacs/gpu_utils/oclutils.h"
62 #include "gromacs/hardware/gpu_hw_info.h"
63 #include "gromacs/math/vectypes.h"
64 #include "gromacs/mdlib/force_flags.h"
65 #include "gromacs/mdtypes/interaction_const.h"
66 #include "gromacs/mdtypes/md_enums.h"
67 #include "gromacs/nbnxm/atomdata.h"
68 #include "gromacs/nbnxm/gpu_data_mgmt.h"
69 #include "gromacs/nbnxm/gpu_jit_support.h"
70 #include "gromacs/nbnxm/nbnxm.h"
71 #include "gromacs/nbnxm/nbnxm_gpu.h"
72 #include "gromacs/nbnxm/pairlistsets.h"
73 #include "gromacs/pbcutil/ishift.h"
74 #include "gromacs/timing/gpu_timing.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/real.h"
79 #include "gromacs/utility/smalloc.h"
81 #include "nbnxm_ocl_internal.h"
83 namespace Nbnxm
86 /*! \brief Copies of values from cl_driver_diagnostics_intel.h,
87 * which isn't guaranteed to be available. */
88 /**@{*/
89 #define CL_CONTEXT_SHOW_DIAGNOSTICS_INTEL 0x4106
90 #define CL_CONTEXT_DIAGNOSTICS_LEVEL_GOOD_INTEL 0x1
91 #define CL_CONTEXT_DIAGNOSTICS_LEVEL_BAD_INTEL 0x2
92 #define CL_CONTEXT_DIAGNOSTICS_LEVEL_NEUTRAL_INTEL 0x4
93 /**@}*/
95 /*! \brief This parameter should be determined heuristically from the
96 * kernel execution times
98 * This value is best for small systems on a single AMD Radeon R9 290X
99 * (and about 5% faster than 40, which is the default for CUDA
100 * devices). Larger simulation systems were quite insensitive to the
101 * value of this parameter.
103 static unsigned int gpu_min_ci_balanced_factor = 50;
106 /*! \brief Returns true if LJ combination rules are used in the non-bonded kernels.
108 * Full doc in nbnxn_ocl_internal.h */
109 bool useLjCombRule(int vdwType)
111 return (vdwType == evdwOclCUTCOMBGEOM || vdwType == evdwOclCUTCOMBLB);
114 /*! \brief Tabulates the Ewald Coulomb force and initializes the size/scale
115 * and the table GPU array.
117 * If called with an already allocated table, it just re-uploads the
118 * table.
120 static void init_ewald_coulomb_force_table(const EwaldCorrectionTables& tables,
121 cl_nbparam_t* nbp,
122 const gmx_device_runtime_data_t* runData)
124 cl_mem coul_tab;
126 cl_int cl_error;
128 if (nbp->coulomb_tab_climg2d != nullptr)
130 freeDeviceBuffer(&(nbp->coulomb_tab_climg2d));
133 /* Switched from using textures to using buffers */
134 // TODO: decide which alternative is most efficient - textures or buffers.
136 cl_image_format array_format;
138 array_format.image_channel_data_type = CL_FLOAT;
139 array_format.image_channel_order = CL_R;
141 coul_tab = clCreateImage2D(runData->context, CL_MEM_READ_WRITE | CL_MEM_COPY_HOST_PTR,
142 &array_format, tabsize, 1, 0, ftmp, &cl_error);
145 coul_tab = clCreateBuffer(
146 runData->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY | CL_MEM_COPY_HOST_PTR,
147 tables.tableF.size() * sizeof(cl_float), const_cast<real*>(tables.tableF.data()), &cl_error);
148 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
149 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
151 nbp->coulomb_tab_climg2d = coul_tab;
152 nbp->coulomb_tab_scale = tables.scale;
156 /*! \brief Initializes the atomdata structure first time, it only gets filled at
157 pair-search.
159 static void init_atomdata_first(cl_atomdata_t* ad, int ntypes, gmx_device_runtime_data_t* runData)
161 cl_int cl_error;
163 ad->ntypes = ntypes;
165 /* An element of the shift_vec device buffer has the same size as one element
166 of the host side shift_vec buffer. */
167 ad->shift_vec_elem_size = sizeof(*nbnxn_atomdata_t::shift_vec.data());
169 ad->shift_vec = clCreateBuffer(runData->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
170 SHIFTS * ad->shift_vec_elem_size, nullptr, &cl_error);
171 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
172 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
173 ad->bShiftVecUploaded = CL_FALSE;
175 /* An element of the fshift device buffer has the same size as one element
176 of the host side fshift buffer. */
177 ad->fshift_elem_size = sizeof(*cl_nb_staging_t::fshift);
179 ad->fshift = clCreateBuffer(runData->context, CL_MEM_READ_WRITE | CL_MEM_HOST_READ_ONLY,
180 SHIFTS * ad->fshift_elem_size, nullptr, &cl_error);
181 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
182 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
184 ad->e_lj = clCreateBuffer(runData->context, CL_MEM_READ_WRITE | CL_MEM_HOST_READ_ONLY,
185 sizeof(float), nullptr, &cl_error);
186 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
187 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
189 ad->e_el = clCreateBuffer(runData->context, CL_MEM_READ_WRITE | CL_MEM_HOST_READ_ONLY,
190 sizeof(float), nullptr, &cl_error);
191 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
192 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
194 /* initialize to nullptr pointers to data that is not allocated here and will
195 need reallocation in nbnxn_gpu_init_atomdata */
196 ad->xq = nullptr;
197 ad->f = nullptr;
199 /* size -1 indicates that the respective array hasn't been initialized yet */
200 ad->natoms = -1;
201 ad->nalloc = -1;
204 /*! \brief Copies all parameters related to the cut-off from ic to nbp
206 static void set_cutoff_parameters(cl_nbparam_t* nbp, const interaction_const_t* ic, const PairlistParams& listParams)
208 nbp->ewald_beta = ic->ewaldcoeff_q;
209 nbp->sh_ewald = ic->sh_ewald;
210 nbp->epsfac = ic->epsfac;
211 nbp->two_k_rf = 2.0 * ic->k_rf;
212 nbp->c_rf = ic->c_rf;
213 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
214 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
215 nbp->rlistOuter_sq = listParams.rlistOuter * listParams.rlistOuter;
216 nbp->rlistInner_sq = listParams.rlistInner * listParams.rlistInner;
217 nbp->useDynamicPruning = listParams.useDynamicPruning;
219 nbp->sh_lj_ewald = ic->sh_lj_ewald;
220 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
222 nbp->rvdw_switch = ic->rvdw_switch;
223 nbp->dispersion_shift = ic->dispersion_shift;
224 nbp->repulsion_shift = ic->repulsion_shift;
225 nbp->vdw_switch = ic->vdw_switch;
228 /*! \brief Returns the kinds of electrostatics and Vdw OpenCL
229 * kernels that will be used.
231 * Respectively, these values are from enum eelOcl and enum
232 * evdwOcl. */
233 static void map_interaction_types_to_gpu_kernel_flavors(const interaction_const_t* ic,
234 int combRule,
235 int* gpu_eeltype,
236 int* gpu_vdwtype)
238 if (ic->vdwtype == evdwCUT)
240 switch (ic->vdw_modifier)
242 case eintmodNONE:
243 case eintmodPOTSHIFT:
244 switch (combRule)
246 case ljcrNONE: *gpu_vdwtype = evdwOclCUT; break;
247 case ljcrGEOM: *gpu_vdwtype = evdwOclCUTCOMBGEOM; break;
248 case ljcrLB: *gpu_vdwtype = evdwOclCUTCOMBLB; break;
249 default:
250 gmx_incons(
251 "The requested LJ combination rule is not implemented in the "
252 "OpenCL GPU accelerated kernels!");
254 break;
255 case eintmodFORCESWITCH: *gpu_vdwtype = evdwOclFSWITCH; break;
256 case eintmodPOTSWITCH: *gpu_vdwtype = evdwOclPSWITCH; break;
257 default:
258 gmx_incons(
259 "The requested VdW interaction modifier is not implemented in the GPU "
260 "accelerated kernels!");
263 else if (ic->vdwtype == evdwPME)
265 if (ic->ljpme_comb_rule == ljcrGEOM)
267 *gpu_vdwtype = evdwOclEWALDGEOM;
269 else
271 *gpu_vdwtype = evdwOclEWALDLB;
274 else
276 gmx_incons("The requested VdW type is not implemented in the GPU accelerated kernels!");
279 if (ic->eeltype == eelCUT)
281 *gpu_eeltype = eelOclCUT;
283 else if (EEL_RF(ic->eeltype))
285 *gpu_eeltype = eelOclRF;
287 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
289 *gpu_eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
291 else
293 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
294 gmx_incons(
295 "The requested electrostatics type is not implemented in the GPU accelerated "
296 "kernels!");
300 /*! \brief Initializes the nonbonded parameter data structure.
302 static void init_nbparam(cl_nbparam_t* nbp,
303 const interaction_const_t* ic,
304 const PairlistParams& listParams,
305 const nbnxn_atomdata_t::Params& nbatParams,
306 const gmx_device_runtime_data_t* runData)
308 cl_int cl_error;
310 set_cutoff_parameters(nbp, ic, listParams);
312 map_interaction_types_to_gpu_kernel_flavors(ic, nbatParams.comb_rule, &(nbp->eeltype), &(nbp->vdwtype));
314 if (ic->vdwtype == evdwPME)
316 if (ic->ljpme_comb_rule == ljcrGEOM)
318 GMX_ASSERT(nbatParams.comb_rule == ljcrGEOM, "Combination rule mismatch!");
320 else
322 GMX_ASSERT(nbatParams.comb_rule == ljcrLB, "Combination rule mismatch!");
325 /* generate table for PME */
326 nbp->coulomb_tab_climg2d = nullptr;
327 if (nbp->eeltype == eelOclEWALD_TAB || nbp->eeltype == eelOclEWALD_TAB_TWIN)
329 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
330 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, runData);
332 else
333 // TODO: improvement needed.
334 // The image2d is created here even if eeltype is not eelCuEWALD_TAB or eelCuEWALD_TAB_TWIN
335 // because the OpenCL kernels don't accept nullptr values for image2D parameters.
337 /* Switched from using textures to using buffers */
338 // TODO: decide which alternative is most efficient - textures or buffers.
340 cl_image_format array_format;
342 array_format.image_channel_data_type = CL_FLOAT;
343 array_format.image_channel_order = CL_R;
345 nbp->coulomb_tab_climg2d = clCreateImage2D(runData->context, CL_MEM_READ_WRITE,
346 &array_format, 1, 1, 0, nullptr, &cl_error);
349 nbp->coulomb_tab_climg2d = clCreateBuffer(runData->context, CL_MEM_READ_ONLY,
350 sizeof(cl_float), nullptr, &cl_error);
351 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
352 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
355 const int nnbfp = 2 * nbatParams.numTypes * nbatParams.numTypes;
356 const int nnbfp_comb = 2 * nbatParams.numTypes;
359 /* Switched from using textures to using buffers */
360 // TODO: decide which alternative is most efficient - textures or buffers.
362 cl_image_format array_format;
364 array_format.image_channel_data_type = CL_FLOAT;
365 array_format.image_channel_order = CL_R;
367 nbp->nbfp_climg2d = clCreateImage2D(runData->context, CL_MEM_READ_ONLY |
368 CL_MEM_COPY_HOST_PTR, &array_format, nnbfp, 1, 0, nbat->nbfp, &cl_error);
371 nbp->nbfp_climg2d = clCreateBuffer(
372 runData->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY | CL_MEM_COPY_HOST_PTR,
373 nnbfp * sizeof(cl_float), const_cast<float*>(nbatParams.nbfp.data()), &cl_error);
374 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
375 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
377 if (ic->vdwtype == evdwPME)
379 /* Switched from using textures to using buffers */
380 // TODO: decide which alternative is most efficient - textures or buffers.
381 /* nbp->nbfp_comb_climg2d = clCreateImage2D(runData->context, CL_MEM_READ_WRITE |
382 CL_MEM_COPY_HOST_PTR, &array_format, nnbfp_comb, 1, 0, nbat->nbfp_comb, &cl_error);*/
383 nbp->nbfp_comb_climg2d = clCreateBuffer(
384 runData->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY | CL_MEM_COPY_HOST_PTR,
385 nnbfp_comb * sizeof(cl_float), const_cast<float*>(nbatParams.nbfp_comb.data()),
386 &cl_error);
387 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
388 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
390 else
392 // TODO: improvement needed.
393 // The image2d is created here even if vdwtype is not evdwPME because the OpenCL kernels
394 // don't accept nullptr values for image2D parameters.
395 /* Switched from using textures to using buffers */
396 // TODO: decide which alternative is most efficient - textures or buffers.
397 /* nbp->nbfp_comb_climg2d = clCreateImage2D(runData->context, CL_MEM_READ_WRITE,
398 &array_format, 1, 1, 0, nullptr, &cl_error);*/
399 nbp->nbfp_comb_climg2d = clCreateBuffer(runData->context, CL_MEM_READ_ONLY,
400 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());
407 //! This function is documented in the header file
408 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t* nbv, const interaction_const_t* ic)
410 if (!nbv || !nbv->useGpu())
412 return;
414 gmx_nbnxn_ocl_t* nb = nbv->gpu_nbv;
415 cl_nbparam_t* nbp = nb->nbparam;
417 set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
419 nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
421 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
422 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, nb->dev_rundata);
425 /*! \brief Initializes the pair list data structure.
427 static void init_plist(cl_plist_t* pl)
429 /* initialize to nullptr pointers to data that is not allocated here and will
430 need reallocation in nbnxn_gpu_init_pairlist */
431 pl->sci = nullptr;
432 pl->cj4 = nullptr;
433 pl->imask = nullptr;
434 pl->excl = nullptr;
436 /* size -1 indicates that the respective array hasn't been initialized yet */
437 pl->na_c = -1;
438 pl->nsci = -1;
439 pl->sci_nalloc = -1;
440 pl->ncj4 = -1;
441 pl->cj4_nalloc = -1;
442 pl->nimask = -1;
443 pl->imask_nalloc = -1;
444 pl->nexcl = -1;
445 pl->excl_nalloc = -1;
446 pl->haveFreshList = false;
449 /*! \brief Initializes the timings data structure.
451 static void init_timings(gmx_wallclock_gpu_nbnxn_t* t)
453 int i, j;
455 t->nb_h2d_t = 0.0;
456 t->nb_d2h_t = 0.0;
457 t->nb_c = 0;
458 t->pl_h2d_t = 0.0;
459 t->pl_h2d_c = 0;
460 for (i = 0; i < 2; i++)
462 for (j = 0; j < 2; j++)
464 t->ktime[i][j].t = 0.0;
465 t->ktime[i][j].c = 0;
469 t->pruneTime.c = 0;
470 t->pruneTime.t = 0.0;
471 t->dynamicPruneTime.c = 0;
472 t->dynamicPruneTime.t = 0.0;
476 //! OpenCL notification callback function
477 static void CL_CALLBACK ocl_notify_fn(const char* pErrInfo,
478 const void gmx_unused* private_info,
479 size_t gmx_unused cb,
480 void gmx_unused* user_data)
482 if (pErrInfo != nullptr)
484 printf("%s\n", pErrInfo); // Print error/hint
488 /*! \brief Creates context for OpenCL GPU given by \p mygpu
490 * A fatal error results if creation fails.
492 * \param[inout] runtimeData runtime data including program and context
493 * \param[in] devInfo device info struct
494 * \param[in] rank MPI rank (for error reporting)
496 static void nbnxn_gpu_create_context(gmx_device_runtime_data_t* runtimeData,
497 const gmx_device_info_t* devInfo,
498 int rank)
500 cl_context_properties context_properties[5];
501 cl_platform_id platform_id;
502 cl_device_id device_id;
503 cl_context context;
504 cl_int cl_error;
506 assert(runtimeData != nullptr);
507 assert(devInfo != nullptr);
509 platform_id = devInfo->ocl_gpu_id.ocl_platform_id;
510 device_id = devInfo->ocl_gpu_id.ocl_device_id;
512 int i = 0;
513 context_properties[i++] = CL_CONTEXT_PLATFORM;
514 context_properties[i++] = reinterpret_cast<cl_context_properties>(platform_id);
515 if (getenv("GMX_OCL_SHOW_DIAGNOSTICS"))
517 context_properties[i++] = CL_CONTEXT_SHOW_DIAGNOSTICS_INTEL;
518 context_properties[i++] =
519 CL_CONTEXT_DIAGNOSTICS_LEVEL_BAD_INTEL | CL_CONTEXT_DIAGNOSTICS_LEVEL_NEUTRAL_INTEL;
521 context_properties[i++] = 0; /* Terminates the list of properties */
523 context = clCreateContext(context_properties, 1, &device_id, ocl_notify_fn, nullptr, &cl_error);
524 if (CL_SUCCESS != cl_error)
526 gmx_fatal(FARGS, "On rank %d failed to create context for GPU #%s:\n OpenCL error %d: %s",
527 rank, devInfo->device_name, cl_error, ocl_get_error_string(cl_error).c_str());
530 runtimeData->context = context;
533 /*! \brief Initializes the OpenCL kernel pointers of the nbnxn_ocl_ptr_t input data structure. */
534 static cl_kernel nbnxn_gpu_create_kernel(gmx_nbnxn_ocl_t* nb, const char* kernel_name)
536 cl_kernel kernel;
537 cl_int cl_error;
539 kernel = clCreateKernel(nb->dev_rundata->program, kernel_name, &cl_error);
540 if (CL_SUCCESS != cl_error)
542 gmx_fatal(FARGS, "Failed to create kernel '%s' for GPU #%s: OpenCL error %d", kernel_name,
543 nb->dev_info->device_name, cl_error);
546 return kernel;
549 /*! \brief Clears nonbonded shift force output array and energy outputs on the GPU.
551 static void nbnxn_ocl_clear_e_fshift(gmx_nbnxn_ocl_t* nb)
554 cl_int cl_error;
555 cl_atomdata_t* adat = nb->atdat;
556 cl_command_queue ls = nb->stream[InteractionLocality::Local];
558 size_t local_work_size[3] = { 1, 1, 1 };
559 size_t global_work_size[3] = { 1, 1, 1 };
561 cl_int shifts = SHIFTS * 3;
563 cl_int arg_no;
565 cl_kernel zero_e_fshift = nb->kernel_zero_e_fshift;
567 local_work_size[0] = 64;
568 // Round the total number of threads up from the array size
569 global_work_size[0] = ((shifts + local_work_size[0] - 1) / local_work_size[0]) * local_work_size[0];
571 arg_no = 0;
572 cl_error = clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_mem), &(adat->fshift));
573 cl_error |= clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_mem), &(adat->e_lj));
574 cl_error |= clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_mem), &(adat->e_el));
575 cl_error |= clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_uint), &shifts);
576 GMX_ASSERT(cl_error == CL_SUCCESS, ocl_get_error_string(cl_error).c_str());
578 cl_error = clEnqueueNDRangeKernel(ls, zero_e_fshift, 3, nullptr, global_work_size,
579 local_work_size, 0, nullptr, nullptr);
580 GMX_ASSERT(cl_error == CL_SUCCESS, ocl_get_error_string(cl_error).c_str());
583 /*! \brief Initializes the OpenCL kernel pointers of the nbnxn_ocl_ptr_t input data structure. */
584 static void nbnxn_gpu_init_kernels(gmx_nbnxn_ocl_t* nb)
586 /* Init to 0 main kernel arrays */
587 /* They will be later on initialized in select_nbnxn_kernel */
588 // TODO: consider always creating all variants of the kernels here so that there is no
589 // need for late call to clCreateKernel -- if that gives any advantage?
590 memset(nb->kernel_ener_noprune_ptr, 0, sizeof(nb->kernel_ener_noprune_ptr));
591 memset(nb->kernel_ener_prune_ptr, 0, sizeof(nb->kernel_ener_prune_ptr));
592 memset(nb->kernel_noener_noprune_ptr, 0, sizeof(nb->kernel_noener_noprune_ptr));
593 memset(nb->kernel_noener_prune_ptr, 0, sizeof(nb->kernel_noener_prune_ptr));
595 /* Init pruning kernels
597 * TODO: we could avoid creating kernels if dynamic pruning is turned off,
598 * but ATM that depends on force flags not passed into the initialization.
600 nb->kernel_pruneonly[epruneFirst] = nbnxn_gpu_create_kernel(nb, "nbnxn_kernel_prune_opencl");
601 nb->kernel_pruneonly[epruneRolling] =
602 nbnxn_gpu_create_kernel(nb, "nbnxn_kernel_prune_rolling_opencl");
604 /* Init auxiliary kernels */
605 nb->kernel_zero_e_fshift = nbnxn_gpu_create_kernel(nb, "zero_e_fshift");
608 /*! \brief Initializes simulation constant data.
610 * Initializes members of the atomdata and nbparam structs and
611 * clears e/fshift output buffers.
613 static void nbnxn_ocl_init_const(gmx_nbnxn_ocl_t* nb,
614 const interaction_const_t* ic,
615 const PairlistParams& listParams,
616 const nbnxn_atomdata_t::Params& nbatParams)
618 init_atomdata_first(nb->atdat, nbatParams.numTypes, nb->dev_rundata);
619 init_nbparam(nb->nbparam, ic, listParams, nbatParams, nb->dev_rundata);
623 //! This function is documented in the header file
624 gmx_nbnxn_ocl_t* gpu_init(const gmx_device_info_t* deviceInfo,
625 const interaction_const_t* ic,
626 const PairlistParams& listParams,
627 const nbnxn_atomdata_t* nbat,
628 const int rank,
629 const gmx_bool bLocalAndNonlocal)
631 gmx_nbnxn_ocl_t* nb;
632 cl_int cl_error;
633 cl_command_queue_properties queue_properties;
635 assert(ic);
637 snew(nb, 1);
638 snew(nb->atdat, 1);
639 snew(nb->nbparam, 1);
640 snew(nb->plist[InteractionLocality::Local], 1);
641 if (bLocalAndNonlocal)
643 snew(nb->plist[InteractionLocality::NonLocal], 1);
646 nb->bUseTwoStreams = static_cast<cl_bool>(bLocalAndNonlocal);
648 nb->timers = new cl_timers_t();
649 snew(nb->timings, 1);
651 /* set device info, just point it to the right GPU among the detected ones */
652 nb->dev_info = deviceInfo;
653 snew(nb->dev_rundata, 1);
655 /* init nbst */
656 pmalloc(reinterpret_cast<void**>(&nb->nbst.e_lj), sizeof(*nb->nbst.e_lj));
657 pmalloc(reinterpret_cast<void**>(&nb->nbst.e_el), sizeof(*nb->nbst.e_el));
658 pmalloc(reinterpret_cast<void**>(&nb->nbst.fshift), SHIFTS * sizeof(*nb->nbst.fshift));
660 init_plist(nb->plist[InteractionLocality::Local]);
662 /* OpenCL timing disabled if GMX_DISABLE_GPU_TIMING is defined. */
663 nb->bDoTime = static_cast<cl_bool>(getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
665 /* Create queues only after bDoTime has been initialized */
666 if (nb->bDoTime)
668 queue_properties = CL_QUEUE_PROFILING_ENABLE;
670 else
672 queue_properties = 0;
675 nbnxn_gpu_create_context(nb->dev_rundata, nb->dev_info, rank);
677 /* local/non-local GPU streams */
678 nb->stream[InteractionLocality::Local] = clCreateCommandQueue(
679 nb->dev_rundata->context, nb->dev_info->ocl_gpu_id.ocl_device_id, queue_properties, &cl_error);
680 if (CL_SUCCESS != cl_error)
682 gmx_fatal(FARGS, "On rank %d failed to create context for GPU #%s: OpenCL error %d", rank,
683 nb->dev_info->device_name, cl_error);
686 if (nb->bUseTwoStreams)
688 init_plist(nb->plist[InteractionLocality::NonLocal]);
690 nb->stream[InteractionLocality::NonLocal] =
691 clCreateCommandQueue(nb->dev_rundata->context, nb->dev_info->ocl_gpu_id.ocl_device_id,
692 queue_properties, &cl_error);
693 if (CL_SUCCESS != cl_error)
695 gmx_fatal(FARGS, "On rank %d failed to create context for GPU #%s: OpenCL error %d",
696 rank, nb->dev_info->device_name, cl_error);
700 if (nb->bDoTime)
702 init_timings(nb->timings);
705 nbnxn_ocl_init_const(nb, ic, listParams, nbat->params());
707 /* Enable LJ param manual prefetch for AMD or Intel or if we request through env. var.
708 * TODO: decide about NVIDIA
710 nb->bPrefetchLjParam = (getenv("GMX_OCL_DISABLE_I_PREFETCH") == nullptr)
711 && ((nb->dev_info->vendor_e == OCL_VENDOR_AMD)
712 || (nb->dev_info->vendor_e == OCL_VENDOR_INTEL)
713 || (getenv("GMX_OCL_ENABLE_I_PREFETCH") != nullptr));
715 /* NOTE: in CUDA we pick L1 cache configuration for the nbnxn kernels here,
716 * but sadly this is not supported in OpenCL (yet?). Consider adding it if
717 * it becomes supported.
719 nbnxn_gpu_compile_kernels(nb);
720 nbnxn_gpu_init_kernels(nb);
722 /* clear energy and shift force outputs */
723 nbnxn_ocl_clear_e_fshift(nb);
725 if (debug)
727 fprintf(debug, "Initialized OpenCL data structures.\n");
730 return nb;
733 /*! \brief Clears the first natoms_clear elements of the GPU nonbonded force output array.
735 static void nbnxn_ocl_clear_f(gmx_nbnxn_ocl_t* nb, int natoms_clear)
737 if (natoms_clear == 0)
739 return;
742 cl_int gmx_used_in_debug cl_error;
744 cl_atomdata_t* atomData = nb->atdat;
745 cl_command_queue ls = nb->stream[InteractionLocality::Local];
746 cl_float value = 0.0F;
748 cl_error = clEnqueueFillBuffer(ls, atomData->f, &value, sizeof(cl_float), 0,
749 natoms_clear * sizeof(rvec), 0, nullptr, nullptr);
750 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
751 ("clEnqueueFillBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
754 //! This function is documented in the header file
755 void gpu_clear_outputs(gmx_nbnxn_ocl_t* nb, bool computeVirial)
757 nbnxn_ocl_clear_f(nb, nb->atdat->natoms);
758 /* clear shift force array and energies if the outputs were
759 used in the current step */
760 if (computeVirial)
762 nbnxn_ocl_clear_e_fshift(nb);
765 /* kick off buffer clearing kernel to ensure concurrency with constraints/update */
766 cl_int gmx_unused cl_error;
767 cl_error = clFlush(nb->stream[InteractionLocality::Local]);
768 assert(CL_SUCCESS == cl_error);
771 //! This function is documented in the header file
772 void gpu_init_pairlist(gmx_nbnxn_ocl_t* nb, const NbnxnPairlistGpu* h_plist, const InteractionLocality iloc)
774 char sbuf[STRLEN];
775 // Timing accumulation should happen only if there was work to do
776 // because getLastRangeTime() gets skipped with empty lists later
777 // which leads to the counter not being reset.
778 bool bDoTime = ((nb->bDoTime == CL_TRUE) && !h_plist->sci.empty());
779 cl_command_queue stream = nb->stream[iloc];
780 cl_plist_t* d_plist = nb->plist[iloc];
782 if (d_plist->na_c < 0)
784 d_plist->na_c = h_plist->na_ci;
786 else
788 if (d_plist->na_c != h_plist->na_ci)
790 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
791 d_plist->na_c, h_plist->na_ci);
792 gmx_incons(sbuf);
796 gpu_timers_t::Interaction& iTimers = nb->timers->interaction[iloc];
798 if (bDoTime)
800 iTimers.pl_h2d.openTimingRegion(stream);
801 iTimers.didPairlistH2D = true;
804 // TODO most of this function is same in CUDA and OpenCL, move into the header
805 DeviceContext context = nb->dev_rundata->context;
807 reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(), &d_plist->nsci, &d_plist->sci_nalloc, context);
808 copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(), stream,
809 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
811 reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(), &d_plist->ncj4, &d_plist->cj4_nalloc, context);
812 copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(), stream,
813 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
815 reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size() * c_nbnxnGpuClusterpairSplit,
816 &d_plist->nimask, &d_plist->imask_nalloc, context);
818 reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(), &d_plist->nexcl,
819 &d_plist->excl_nalloc, context);
820 copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(), stream,
821 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
823 if (bDoTime)
825 iTimers.pl_h2d.closeTimingRegion(stream);
828 /* need to prune the pair list during the next step */
829 d_plist->haveFreshList = true;
832 //! This function is documented in the header file
833 void gpu_upload_shiftvec(gmx_nbnxn_ocl_t* nb, const nbnxn_atomdata_t* nbatom)
835 cl_atomdata_t* adat = nb->atdat;
836 cl_command_queue ls = nb->stream[InteractionLocality::Local];
838 /* only if we have a dynamic box */
839 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
841 ocl_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(), 0,
842 SHIFTS * adat->shift_vec_elem_size, ls, nullptr);
843 adat->bShiftVecUploaded = CL_TRUE;
847 //! This function is documented in the header file
848 void gpu_init_atomdata(gmx_nbnxn_ocl_t* nb, const nbnxn_atomdata_t* nbat)
850 cl_int cl_error;
851 int nalloc, natoms;
852 bool realloced;
853 bool bDoTime = nb->bDoTime == CL_TRUE;
854 cl_timers_t* timers = nb->timers;
855 cl_atomdata_t* d_atdat = nb->atdat;
856 cl_command_queue ls = nb->stream[InteractionLocality::Local];
858 natoms = nbat->numAtoms();
859 realloced = false;
861 if (bDoTime)
863 /* time async copy */
864 timers->atdat.openTimingRegion(ls);
867 /* need to reallocate if we have to copy more atoms than the amount of space
868 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
869 if (natoms > d_atdat->nalloc)
871 nalloc = over_alloc_small(natoms);
873 /* free up first if the arrays have already been initialized */
874 if (d_atdat->nalloc != -1)
876 freeDeviceBuffer(&d_atdat->f);
877 freeDeviceBuffer(&d_atdat->xq);
878 freeDeviceBuffer(&d_atdat->lj_comb);
879 freeDeviceBuffer(&d_atdat->atom_types);
882 d_atdat->f_elem_size = sizeof(rvec);
884 d_atdat->f = clCreateBuffer(nb->dev_rundata->context, CL_MEM_READ_WRITE | CL_MEM_HOST_READ_ONLY,
885 nalloc * d_atdat->f_elem_size, nullptr, &cl_error);
886 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
887 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
889 d_atdat->xq = clCreateBuffer(nb->dev_rundata->context, CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
890 nalloc * sizeof(cl_float4), nullptr, &cl_error);
891 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
892 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
894 if (useLjCombRule(nb->nbparam->vdwtype))
896 d_atdat->lj_comb = clCreateBuffer(nb->dev_rundata->context,
897 CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
898 nalloc * sizeof(cl_float2), nullptr, &cl_error);
899 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
900 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
902 else
904 d_atdat->atom_types = clCreateBuffer(nb->dev_rundata->context,
905 CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
906 nalloc * sizeof(int), nullptr, &cl_error);
907 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
908 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
911 d_atdat->nalloc = nalloc;
912 realloced = true;
915 d_atdat->natoms = natoms;
916 d_atdat->natoms_local = nbat->natoms_local;
918 /* need to clear GPU f output if realloc happened */
919 if (realloced)
921 nbnxn_ocl_clear_f(nb, nalloc);
924 if (useLjCombRule(nb->nbparam->vdwtype))
926 ocl_copy_H2D_async(d_atdat->lj_comb, nbat->params().lj_comb.data(), 0, natoms * sizeof(cl_float2),
927 ls, bDoTime ? timers->atdat.fetchNextEvent() : nullptr);
929 else
931 ocl_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(), 0, natoms * sizeof(int),
932 ls, bDoTime ? timers->atdat.fetchNextEvent() : nullptr);
935 if (bDoTime)
937 timers->atdat.closeTimingRegion(ls);
940 /* kick off the tasks enqueued above to ensure concurrency with the search */
941 cl_error = clFlush(ls);
942 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
943 ("clFlush failed: " + ocl_get_error_string(cl_error)).c_str());
946 /*! \brief Releases an OpenCL kernel pointer */
947 static void free_kernel(cl_kernel* kernel_ptr)
949 cl_int gmx_unused cl_error;
951 assert(nullptr != kernel_ptr);
953 if (*kernel_ptr)
955 cl_error = clReleaseKernel(*kernel_ptr);
956 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
957 ("clReleaseKernel failed: " + ocl_get_error_string(cl_error)).c_str());
959 *kernel_ptr = nullptr;
963 /*! \brief Releases a list of OpenCL kernel pointers */
964 static void free_kernels(cl_kernel* kernels, int count)
966 int i;
968 for (i = 0; i < count; i++)
970 free_kernel(kernels + i);
974 /*! \brief Free the OpenCL runtime data (context and program).
976 * The function releases the OpenCL context and program assuciated with the
977 * device that the calling PP rank is running on.
979 * \param runData [in] porinter to the structure with runtime data.
981 static void free_gpu_device_runtime_data(gmx_device_runtime_data_t* runData)
983 if (runData == nullptr)
985 return;
988 cl_int gmx_unused cl_error;
990 if (runData->context)
992 cl_error = clReleaseContext(runData->context);
993 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
994 ("clReleaseContext failed: " + ocl_get_error_string(cl_error)).c_str());
995 runData->context = nullptr;
998 if (runData->program)
1000 cl_error = clReleaseProgram(runData->program);
1001 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
1002 ("clReleaseProgram failed: " + ocl_get_error_string(cl_error)).c_str());
1003 runData->program = nullptr;
1007 //! This function is documented in the header file
1008 void gpu_free(gmx_nbnxn_ocl_t* nb)
1010 if (nb == nullptr)
1012 return;
1015 /* Free kernels */
1016 int kernel_count = sizeof(nb->kernel_ener_noprune_ptr) / sizeof(nb->kernel_ener_noprune_ptr[0][0]);
1017 free_kernels(nb->kernel_ener_noprune_ptr[0], kernel_count);
1019 kernel_count = sizeof(nb->kernel_ener_prune_ptr) / sizeof(nb->kernel_ener_prune_ptr[0][0]);
1020 free_kernels(nb->kernel_ener_prune_ptr[0], kernel_count);
1022 kernel_count = sizeof(nb->kernel_noener_noprune_ptr) / sizeof(nb->kernel_noener_noprune_ptr[0][0]);
1023 free_kernels(nb->kernel_noener_noprune_ptr[0], kernel_count);
1025 kernel_count = sizeof(nb->kernel_noener_prune_ptr) / sizeof(nb->kernel_noener_prune_ptr[0][0]);
1026 free_kernels(nb->kernel_noener_prune_ptr[0], kernel_count);
1028 free_kernel(&(nb->kernel_zero_e_fshift));
1030 /* Free atdat */
1031 freeDeviceBuffer(&(nb->atdat->xq));
1032 freeDeviceBuffer(&(nb->atdat->f));
1033 freeDeviceBuffer(&(nb->atdat->e_lj));
1034 freeDeviceBuffer(&(nb->atdat->e_el));
1035 freeDeviceBuffer(&(nb->atdat->fshift));
1036 freeDeviceBuffer(&(nb->atdat->lj_comb));
1037 freeDeviceBuffer(&(nb->atdat->atom_types));
1038 freeDeviceBuffer(&(nb->atdat->shift_vec));
1039 sfree(nb->atdat);
1041 /* Free nbparam */
1042 freeDeviceBuffer(&(nb->nbparam->nbfp_climg2d));
1043 freeDeviceBuffer(&(nb->nbparam->nbfp_comb_climg2d));
1044 freeDeviceBuffer(&(nb->nbparam->coulomb_tab_climg2d));
1045 sfree(nb->nbparam);
1047 /* Free plist */
1048 auto* plist = nb->plist[InteractionLocality::Local];
1049 freeDeviceBuffer(&plist->sci);
1050 freeDeviceBuffer(&plist->cj4);
1051 freeDeviceBuffer(&plist->imask);
1052 freeDeviceBuffer(&plist->excl);
1053 sfree(plist);
1054 if (nb->bUseTwoStreams)
1056 auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
1057 freeDeviceBuffer(&plist_nl->sci);
1058 freeDeviceBuffer(&plist_nl->cj4);
1059 freeDeviceBuffer(&plist_nl->imask);
1060 freeDeviceBuffer(&plist_nl->excl);
1061 sfree(plist_nl);
1064 /* Free nbst */
1065 pfree(nb->nbst.e_lj);
1066 nb->nbst.e_lj = nullptr;
1068 pfree(nb->nbst.e_el);
1069 nb->nbst.e_el = nullptr;
1071 pfree(nb->nbst.fshift);
1072 nb->nbst.fshift = nullptr;
1074 /* Free command queues */
1075 clReleaseCommandQueue(nb->stream[InteractionLocality::Local]);
1076 nb->stream[InteractionLocality::Local] = nullptr;
1077 if (nb->bUseTwoStreams)
1079 clReleaseCommandQueue(nb->stream[InteractionLocality::NonLocal]);
1080 nb->stream[InteractionLocality::NonLocal] = nullptr;
1082 /* Free other events */
1083 if (nb->nonlocal_done)
1085 clReleaseEvent(nb->nonlocal_done);
1086 nb->nonlocal_done = nullptr;
1088 if (nb->misc_ops_and_local_H2D_done)
1090 clReleaseEvent(nb->misc_ops_and_local_H2D_done);
1091 nb->misc_ops_and_local_H2D_done = nullptr;
1094 free_gpu_device_runtime_data(nb->dev_rundata);
1095 sfree(nb->dev_rundata);
1097 /* Free timers and timings */
1098 delete nb->timers;
1099 sfree(nb->timings);
1100 sfree(nb);
1102 if (debug)
1104 fprintf(debug, "Cleaned up OpenCL data structures.\n");
1108 //! This function is documented in the header file
1109 gmx_wallclock_gpu_nbnxn_t* gpu_get_timings(gmx_nbnxn_ocl_t* nb)
1111 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
1114 //! This function is documented in the header file
1115 void gpu_reset_timings(nonbonded_verlet_t* nbv)
1117 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
1119 init_timings(nbv->gpu_nbv->timings);
1123 //! This function is documented in the header file
1124 int gpu_min_ci_balanced(gmx_nbnxn_ocl_t* nb)
1126 return nb != nullptr ? gpu_min_ci_balanced_factor * nb->dev_info->compute_units : 0;
1129 //! This function is documented in the header file
1130 gmx_bool gpu_is_kernel_ewald_analytical(const gmx_nbnxn_ocl_t* nb)
1132 return ((nb->nbparam->eeltype == eelOclEWALD_ANA) || (nb->nbparam->eeltype == eelOclEWALD_ANA_TWIN));
1135 } // namespace Nbnxm