Use init(..) function to build DeviceContext
[gromacs.git] / src / gromacs / nbnxm / opencl / nbnxm_ocl_data_mgmt.cpp
blob5b5911941d86f1d7bd8467b5f823a81f9f35e401
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 #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"
76 #include "nbnxm_ocl_internal.h"
77 #include "nbnxm_ocl_types.h"
79 namespace Nbnxm
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 /**@}*/
91 /*! \brief This parameter should be determined heuristically from the
92 * kernel execution times
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.
99 static unsigned int gpu_min_ci_balanced_factor = 50;
102 /*! \brief Returns true if LJ combination rules are used in the non-bonded kernels.
104 * Full doc in nbnxn_ocl_internal.h */
105 bool useLjCombRule(int vdwType)
107 return (vdwType == evdwOclCUTCOMBGEOM || vdwType == evdwOclCUTCOMBLB);
110 /*! \brief Tabulates the Ewald Coulomb force and initializes the size/scale
111 * and the table GPU array.
113 * If called with an already allocated table, it just re-uploads the
114 * table.
116 static void init_ewald_coulomb_force_table(const EwaldCorrectionTables& tables,
117 cl_nbparam_t* nbp,
118 const gmx_device_runtime_data_t* runData)
120 cl_mem coul_tab;
122 cl_int cl_error;
124 if (nbp->coulomb_tab_climg2d != nullptr)
126 freeDeviceBuffer(&(nbp->coulomb_tab_climg2d));
129 /* Switched from using textures to using buffers */
130 // TODO: decide which alternative is most efficient - textures or buffers.
132 cl_image_format array_format;
134 array_format.image_channel_data_type = CL_FLOAT;
135 array_format.image_channel_order = CL_R;
137 coul_tab = clCreateImage2D(runData->deviceContext.context(), CL_MEM_READ_WRITE |
138 CL_MEM_COPY_HOST_PTR, &array_format, tabsize, 1, 0, ftmp, &cl_error);
141 coul_tab = clCreateBuffer(runData->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());
148 nbp->coulomb_tab_climg2d = coul_tab;
149 nbp->coulomb_tab_scale = tables.scale;
153 /*! \brief Initializes the atomdata structure first time, it only gets filled at
154 pair-search.
156 static void init_atomdata_first(cl_atomdata_t* ad, int ntypes, gmx_device_runtime_data_t* runData)
158 cl_int cl_error;
160 ad->ntypes = ntypes;
162 ad->shift_vec =
163 clCreateBuffer(runData->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;
169 ad->fshift = clCreateBuffer(runData->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());
174 ad->e_lj = clCreateBuffer(runData->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());
179 ad->e_el = clCreateBuffer(runData->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());
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;
189 /* size -1 indicates that the respective array hasn't been initialized yet */
190 ad->natoms = -1;
191 ad->nalloc = -1;
194 /*! \brief Copies all parameters related to the cut-off from ic to nbp
196 static void set_cutoff_parameters(cl_nbparam_t* nbp, const interaction_const_t* ic, const PairlistParams& listParams)
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;
209 nbp->sh_lj_ewald = ic->sh_lj_ewald;
210 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
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;
218 /*! \brief Returns the kinds of electrostatics and Vdw OpenCL
219 * kernels that will be used.
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)
228 if (ic->vdwtype == evdwCUT)
230 switch (ic->vdw_modifier)
232 case eintmodNONE:
233 case eintmodPOTSHIFT:
234 switch (combRule)
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!");
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!");
253 else if (ic->vdwtype == evdwPME)
255 if (ic->ljpme_comb_rule == ljcrGEOM)
257 *gpu_vdwtype = evdwOclEWALDGEOM;
259 else
261 *gpu_vdwtype = evdwOclEWALDLB;
264 else
266 gmx_incons("The requested VdW type is not implemented in the GPU accelerated kernels!");
269 if (ic->eeltype == eelCUT)
271 *gpu_eeltype = eelOclCUT;
273 else if (EEL_RF(ic->eeltype))
275 *gpu_eeltype = eelOclRF;
277 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
279 *gpu_eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
281 else
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!");
290 /*! \brief Initializes the nonbonded parameter data structure.
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 gmx_device_runtime_data_t* runData)
298 cl_int cl_error;
300 set_cutoff_parameters(nbp, ic, listParams);
302 map_interaction_types_to_gpu_kernel_flavors(ic, nbatParams.comb_rule, &(nbp->eeltype), &(nbp->vdwtype));
304 if (ic->vdwtype == evdwPME)
306 if (ic->ljpme_comb_rule == ljcrGEOM)
308 GMX_ASSERT(nbatParams.comb_rule == ljcrGEOM, "Combination rule mismatch!");
310 else
312 GMX_ASSERT(nbatParams.comb_rule == ljcrLB, "Combination rule mismatch!");
315 /* generate table for PME */
316 nbp->coulomb_tab_climg2d = nullptr;
317 if (nbp->eeltype == eelOclEWALD_TAB || nbp->eeltype == eelOclEWALD_TAB_TWIN)
319 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
320 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, runData);
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.
327 /* Switched from using textures to using buffers */
328 // TODO: decide which alternative is most efficient - textures or buffers.
330 cl_image_format array_format;
332 array_format.image_channel_data_type = CL_FLOAT;
333 array_format.image_channel_order = CL_R;
335 nbp->coulomb_tab_climg2d = clCreateImage2D(runData->deviceContext.context(),
336 CL_MEM_READ_WRITE, &array_format, 1, 1, 0, nullptr, &cl_error);
339 nbp->coulomb_tab_climg2d = clCreateBuffer(runData->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());
345 const int nnbfp = 2 * nbatParams.numTypes * nbatParams.numTypes;
346 const int nnbfp_comb = 2 * nbatParams.numTypes;
349 /* Switched from using textures to using buffers */
350 // TODO: decide which alternative is most efficient - textures or buffers.
352 cl_image_format array_format;
354 array_format.image_channel_data_type = CL_FLOAT;
355 array_format.image_channel_order = CL_R;
357 nbp->nbfp_climg2d = clCreateImage2D(runData->deviceContext.context(), CL_MEM_READ_ONLY |
358 CL_MEM_COPY_HOST_PTR, &array_format, nnbfp, 1, 0, nbat->nbfp, &cl_error);
361 nbp->nbfp_climg2d = clCreateBuffer(
362 runData->deviceContext.context(),
363 CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY | CL_MEM_COPY_HOST_PTR,
364 nnbfp * sizeof(cl_float), const_cast<float*>(nbatParams.nbfp.data()), &cl_error);
365 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
366 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
368 if (ic->vdwtype == evdwPME)
370 /* Switched from using textures to using buffers */
371 // TODO: decide which alternative is most efficient - textures or buffers.
372 /* nbp->nbfp_comb_climg2d = clCreateImage2D(runData->deviceContext.context(), CL_MEM_READ_WRITE |
373 CL_MEM_COPY_HOST_PTR, &array_format, nnbfp_comb, 1, 0, nbat->nbfp_comb, &cl_error);*/
374 nbp->nbfp_comb_climg2d =
375 clCreateBuffer(runData->deviceContext.context(),
376 CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY | CL_MEM_COPY_HOST_PTR,
377 nnbfp_comb * sizeof(cl_float),
378 const_cast<float*>(nbatParams.nbfp_comb.data()), &cl_error);
379 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
380 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
382 else
384 // TODO: improvement needed.
385 // The image2d is created here even if vdwtype is not evdwPME because the OpenCL kernels
386 // don't accept nullptr values for image2D parameters.
387 /* Switched from using textures to using buffers */
388 // TODO: decide which alternative is most efficient - textures or buffers.
389 /* nbp->nbfp_comb_climg2d = clCreateImage2D(runData->deviceContext.context(),
390 CL_MEM_READ_WRITE, &array_format, 1, 1, 0, nullptr, &cl_error);*/
391 nbp->nbfp_comb_climg2d = clCreateBuffer(runData->deviceContext.context(), CL_MEM_READ_ONLY,
392 sizeof(cl_float), nullptr, &cl_error);
393 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
394 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
399 //! This function is documented in the header file
400 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t* nbv, const interaction_const_t* ic)
402 if (!nbv || !nbv->useGpu())
404 return;
406 NbnxmGpu* nb = nbv->gpu_nbv;
407 cl_nbparam_t* nbp = nb->nbparam;
409 set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
411 nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
413 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
414 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, nb->dev_rundata);
417 /*! \brief Initializes the pair list data structure.
419 static void init_plist(cl_plist_t* pl)
421 /* initialize to nullptr pointers to data that is not allocated here and will
422 need reallocation in nbnxn_gpu_init_pairlist */
423 pl->sci = nullptr;
424 pl->cj4 = nullptr;
425 pl->imask = nullptr;
426 pl->excl = nullptr;
428 /* size -1 indicates that the respective array hasn't been initialized yet */
429 pl->na_c = -1;
430 pl->nsci = -1;
431 pl->sci_nalloc = -1;
432 pl->ncj4 = -1;
433 pl->cj4_nalloc = -1;
434 pl->nimask = -1;
435 pl->imask_nalloc = -1;
436 pl->nexcl = -1;
437 pl->excl_nalloc = -1;
438 pl->haveFreshList = false;
441 /*! \brief Initializes the timings data structure.
443 static void init_timings(gmx_wallclock_gpu_nbnxn_t* t)
445 int i, j;
447 t->nb_h2d_t = 0.0;
448 t->nb_d2h_t = 0.0;
449 t->nb_c = 0;
450 t->pl_h2d_t = 0.0;
451 t->pl_h2d_c = 0;
452 for (i = 0; i < 2; i++)
454 for (j = 0; j < 2; j++)
456 t->ktime[i][j].t = 0.0;
457 t->ktime[i][j].c = 0;
461 t->pruneTime.c = 0;
462 t->pruneTime.t = 0.0;
463 t->dynamicPruneTime.c = 0;
464 t->dynamicPruneTime.t = 0.0;
467 /*! \brief Initializes the OpenCL kernel pointers of the nbnxn_ocl_ptr_t input data structure. */
468 static cl_kernel nbnxn_gpu_create_kernel(NbnxmGpu* nb, const char* kernel_name)
470 cl_kernel kernel;
471 cl_int cl_error;
473 kernel = clCreateKernel(nb->dev_rundata->program, kernel_name, &cl_error);
474 if (CL_SUCCESS != cl_error)
476 gmx_fatal(FARGS, "Failed to create kernel '%s' for GPU #%s: OpenCL error %d", kernel_name,
477 nb->deviceInfo->device_name, cl_error);
480 return kernel;
483 /*! \brief Clears nonbonded shift force output array and energy outputs on the GPU.
485 static void nbnxn_ocl_clear_e_fshift(NbnxmGpu* nb)
488 cl_int cl_error;
489 cl_atomdata_t* adat = nb->atdat;
490 cl_command_queue ls = nb->stream[InteractionLocality::Local];
492 size_t local_work_size[3] = { 1, 1, 1 };
493 size_t global_work_size[3] = { 1, 1, 1 };
495 cl_int shifts = SHIFTS * 3;
497 cl_int arg_no;
499 cl_kernel zero_e_fshift = nb->kernel_zero_e_fshift;
501 local_work_size[0] = 64;
502 // Round the total number of threads up from the array size
503 global_work_size[0] = ((shifts + local_work_size[0] - 1) / local_work_size[0]) * local_work_size[0];
505 arg_no = 0;
506 cl_error = clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_mem), &(adat->fshift));
507 cl_error |= clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_mem), &(adat->e_lj));
508 cl_error |= clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_mem), &(adat->e_el));
509 cl_error |= clSetKernelArg(zero_e_fshift, arg_no++, sizeof(cl_uint), &shifts);
510 GMX_ASSERT(cl_error == CL_SUCCESS, ocl_get_error_string(cl_error).c_str());
512 cl_error = clEnqueueNDRangeKernel(ls, zero_e_fshift, 3, nullptr, global_work_size,
513 local_work_size, 0, nullptr, nullptr);
514 GMX_ASSERT(cl_error == CL_SUCCESS, ocl_get_error_string(cl_error).c_str());
517 /*! \brief Initializes the OpenCL kernel pointers of the nbnxn_ocl_ptr_t input data structure. */
518 static void nbnxn_gpu_init_kernels(NbnxmGpu* nb)
520 /* Init to 0 main kernel arrays */
521 /* They will be later on initialized in select_nbnxn_kernel */
522 // TODO: consider always creating all variants of the kernels here so that there is no
523 // need for late call to clCreateKernel -- if that gives any advantage?
524 memset(nb->kernel_ener_noprune_ptr, 0, sizeof(nb->kernel_ener_noprune_ptr));
525 memset(nb->kernel_ener_prune_ptr, 0, sizeof(nb->kernel_ener_prune_ptr));
526 memset(nb->kernel_noener_noprune_ptr, 0, sizeof(nb->kernel_noener_noprune_ptr));
527 memset(nb->kernel_noener_prune_ptr, 0, sizeof(nb->kernel_noener_prune_ptr));
529 /* Init pruning kernels
531 * TODO: we could avoid creating kernels if dynamic pruning is turned off,
532 * but ATM that depends on force flags not passed into the initialization.
534 nb->kernel_pruneonly[epruneFirst] = nbnxn_gpu_create_kernel(nb, "nbnxn_kernel_prune_opencl");
535 nb->kernel_pruneonly[epruneRolling] =
536 nbnxn_gpu_create_kernel(nb, "nbnxn_kernel_prune_rolling_opencl");
538 /* Init auxiliary kernels */
539 nb->kernel_zero_e_fshift = nbnxn_gpu_create_kernel(nb, "zero_e_fshift");
542 /*! \brief Initializes simulation constant data.
544 * Initializes members of the atomdata and nbparam structs and
545 * clears e/fshift output buffers.
547 static void nbnxn_ocl_init_const(NbnxmGpu* nb,
548 const interaction_const_t* ic,
549 const PairlistParams& listParams,
550 const nbnxn_atomdata_t::Params& nbatParams)
552 init_atomdata_first(nb->atdat, nbatParams.numTypes, nb->dev_rundata);
553 init_nbparam(nb->nbparam, ic, listParams, nbatParams, nb->dev_rundata);
557 //! This function is documented in the header file
558 NbnxmGpu* gpu_init(const DeviceInformation* deviceInfo,
559 const interaction_const_t* ic,
560 const PairlistParams& listParams,
561 const nbnxn_atomdata_t* nbat,
562 const int rank,
563 const bool bLocalAndNonlocal)
565 cl_int cl_error;
566 cl_command_queue_properties queue_properties;
568 GMX_ASSERT(ic, "Need a valid interaction constants object");
570 auto nb = new NbnxmGpu;
571 snew(nb->atdat, 1);
572 snew(nb->nbparam, 1);
573 snew(nb->plist[InteractionLocality::Local], 1);
574 if (bLocalAndNonlocal)
576 snew(nb->plist[InteractionLocality::NonLocal], 1);
579 nb->bUseTwoStreams = bLocalAndNonlocal;
581 nb->timers = new cl_timers_t();
582 snew(nb->timings, 1);
584 /* set device info, just point it to the right GPU among the detected ones */
585 nb->deviceInfo = deviceInfo;
586 nb->dev_rundata = new gmx_device_runtime_data_t();
588 /* init nbst */
589 pmalloc(reinterpret_cast<void**>(&nb->nbst.e_lj), sizeof(*nb->nbst.e_lj));
590 pmalloc(reinterpret_cast<void**>(&nb->nbst.e_el), sizeof(*nb->nbst.e_el));
591 pmalloc(reinterpret_cast<void**>(&nb->nbst.fshift), SHIFTS * sizeof(*nb->nbst.fshift));
593 init_plist(nb->plist[InteractionLocality::Local]);
595 /* OpenCL timing disabled if GMX_DISABLE_GPU_TIMING is defined. */
596 nb->bDoTime = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
598 /* Create queues only after bDoTime has been initialized */
599 if (nb->bDoTime)
601 queue_properties = CL_QUEUE_PROFILING_ENABLE;
603 else
605 queue_properties = 0;
608 nb->dev_rundata->deviceContext.init(*deviceInfo);
610 /* local/non-local GPU streams */
611 nb->stream[InteractionLocality::Local] =
612 clCreateCommandQueue(nb->dev_rundata->deviceContext.context(),
613 nb->deviceInfo->oclDeviceId, queue_properties, &cl_error);
614 if (CL_SUCCESS != cl_error)
616 gmx_fatal(FARGS, "On rank %d failed to create context for GPU #%s: OpenCL error %d", rank,
617 nb->deviceInfo->device_name, cl_error);
620 if (nb->bUseTwoStreams)
622 init_plist(nb->plist[InteractionLocality::NonLocal]);
624 nb->stream[InteractionLocality::NonLocal] =
625 clCreateCommandQueue(nb->dev_rundata->deviceContext.context(),
626 nb->deviceInfo->oclDeviceId, queue_properties, &cl_error);
627 if (CL_SUCCESS != cl_error)
629 gmx_fatal(FARGS, "On rank %d failed to create context for GPU #%s: OpenCL error %d",
630 rank, nb->deviceInfo->device_name, cl_error);
634 if (nb->bDoTime)
636 init_timings(nb->timings);
639 nbnxn_ocl_init_const(nb, ic, listParams, nbat->params());
641 /* Enable LJ param manual prefetch for AMD or Intel or if we request through env. var.
642 * TODO: decide about NVIDIA
644 nb->bPrefetchLjParam = (getenv("GMX_OCL_DISABLE_I_PREFETCH") == nullptr)
645 && ((nb->deviceInfo->deviceVendor == DeviceVendor::Amd)
646 || (nb->deviceInfo->deviceVendor == DeviceVendor::Intel)
647 || (getenv("GMX_OCL_ENABLE_I_PREFETCH") != nullptr));
649 /* NOTE: in CUDA we pick L1 cache configuration for the nbnxn kernels here,
650 * but sadly this is not supported in OpenCL (yet?). Consider adding it if
651 * it becomes supported.
653 nbnxn_gpu_compile_kernels(nb);
654 nbnxn_gpu_init_kernels(nb);
656 /* clear energy and shift force outputs */
657 nbnxn_ocl_clear_e_fshift(nb);
659 if (debug)
661 fprintf(debug, "Initialized OpenCL data structures.\n");
664 return nb;
667 /*! \brief Clears the first natoms_clear elements of the GPU nonbonded force output array.
669 static void nbnxn_ocl_clear_f(NbnxmGpu* nb, int natoms_clear)
671 if (natoms_clear == 0)
673 return;
676 cl_int gmx_used_in_debug cl_error;
678 cl_atomdata_t* atomData = nb->atdat;
679 cl_command_queue ls = nb->stream[InteractionLocality::Local];
680 cl_float value = 0.0F;
682 cl_error = clEnqueueFillBuffer(ls, atomData->f, &value, sizeof(cl_float), 0,
683 natoms_clear * sizeof(rvec), 0, nullptr, nullptr);
684 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
685 ("clEnqueueFillBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
688 //! This function is documented in the header file
689 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
691 nbnxn_ocl_clear_f(nb, nb->atdat->natoms);
692 /* clear shift force array and energies if the outputs were
693 used in the current step */
694 if (computeVirial)
696 nbnxn_ocl_clear_e_fshift(nb);
699 /* kick off buffer clearing kernel to ensure concurrency with constraints/update */
700 cl_int gmx_unused cl_error;
701 cl_error = clFlush(nb->stream[InteractionLocality::Local]);
702 GMX_ASSERT(cl_error == CL_SUCCESS, ("clFlush failed: " + ocl_get_error_string(cl_error)).c_str());
705 //! This function is documented in the header file
706 void gpu_init_pairlist(NbnxmGpu* nb, const NbnxnPairlistGpu* h_plist, const InteractionLocality iloc)
708 char sbuf[STRLEN];
709 // Timing accumulation should happen only if there was work to do
710 // because getLastRangeTime() gets skipped with empty lists later
711 // which leads to the counter not being reset.
712 bool bDoTime = (nb->bDoTime && !h_plist->sci.empty());
713 cl_command_queue stream = nb->stream[iloc];
714 cl_plist_t* d_plist = nb->plist[iloc];
716 if (d_plist->na_c < 0)
718 d_plist->na_c = h_plist->na_ci;
720 else
722 if (d_plist->na_c != h_plist->na_ci)
724 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
725 d_plist->na_c, h_plist->na_ci);
726 gmx_incons(sbuf);
730 gpu_timers_t::Interaction& iTimers = nb->timers->interaction[iloc];
732 if (bDoTime)
734 iTimers.pl_h2d.openTimingRegion(stream);
735 iTimers.didPairlistH2D = true;
738 // TODO most of this function is same in CUDA and OpenCL, move into the header
739 const DeviceContext& deviceContext = nb->dev_rundata->deviceContext;
741 reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(), &d_plist->nsci, &d_plist->sci_nalloc,
742 deviceContext);
743 copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(), stream,
744 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
746 reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(), &d_plist->ncj4, &d_plist->cj4_nalloc,
747 deviceContext);
748 copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(), stream,
749 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
751 reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size() * c_nbnxnGpuClusterpairSplit,
752 &d_plist->nimask, &d_plist->imask_nalloc, deviceContext);
754 reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(), &d_plist->nexcl,
755 &d_plist->excl_nalloc, deviceContext);
756 copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(), stream,
757 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
759 if (bDoTime)
761 iTimers.pl_h2d.closeTimingRegion(stream);
764 /* need to prune the pair list during the next step */
765 d_plist->haveFreshList = true;
768 //! This function is documented in the header file
769 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
771 cl_atomdata_t* adat = nb->atdat;
772 cl_command_queue ls = nb->stream[InteractionLocality::Local];
774 /* only if we have a dynamic box */
775 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
777 ocl_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(), 0,
778 SHIFTS * sizeof(nbatom->shift_vec[0]), ls, nullptr);
779 adat->bShiftVecUploaded = CL_TRUE;
783 //! This function is documented in the header file
784 void gpu_init_atomdata(NbnxmGpu* nb, const nbnxn_atomdata_t* nbat)
786 cl_int cl_error;
787 int nalloc, natoms;
788 bool realloced;
789 bool bDoTime = nb->bDoTime;
790 cl_timers_t* timers = nb->timers;
791 cl_atomdata_t* d_atdat = nb->atdat;
792 cl_command_queue ls = nb->stream[InteractionLocality::Local];
794 natoms = nbat->numAtoms();
795 realloced = false;
797 if (bDoTime)
799 /* time async copy */
800 timers->atdat.openTimingRegion(ls);
803 /* need to reallocate if we have to copy more atoms than the amount of space
804 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
805 if (natoms > d_atdat->nalloc)
807 nalloc = over_alloc_small(natoms);
809 /* free up first if the arrays have already been initialized */
810 if (d_atdat->nalloc != -1)
812 freeDeviceBuffer(&d_atdat->f);
813 freeDeviceBuffer(&d_atdat->xq);
814 freeDeviceBuffer(&d_atdat->lj_comb);
815 freeDeviceBuffer(&d_atdat->atom_types);
818 d_atdat->f = clCreateBuffer(nb->dev_rundata->deviceContext.context(),
819 CL_MEM_READ_WRITE | CL_MEM_HOST_READ_ONLY,
820 nalloc * DIM * sizeof(nbat->out[0].f[0]), nullptr, &cl_error);
821 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
822 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
824 d_atdat->xq = clCreateBuffer(nb->dev_rundata->deviceContext.context(),
825 CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
826 nalloc * sizeof(cl_float4), nullptr, &cl_error);
827 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
828 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
830 if (useLjCombRule(nb->nbparam->vdwtype))
832 d_atdat->lj_comb = clCreateBuffer(nb->dev_rundata->deviceContext.context(),
833 CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
834 nalloc * sizeof(cl_float2), nullptr, &cl_error);
835 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
836 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
838 else
840 d_atdat->atom_types = clCreateBuffer(nb->dev_rundata->deviceContext.context(),
841 CL_MEM_READ_ONLY | CL_MEM_HOST_WRITE_ONLY,
842 nalloc * sizeof(int), nullptr, &cl_error);
843 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
844 ("clCreateBuffer failed: " + ocl_get_error_string(cl_error)).c_str());
847 d_atdat->nalloc = nalloc;
848 realloced = true;
851 d_atdat->natoms = natoms;
852 d_atdat->natoms_local = nbat->natoms_local;
854 /* need to clear GPU f output if realloc happened */
855 if (realloced)
857 nbnxn_ocl_clear_f(nb, nalloc);
860 if (useLjCombRule(nb->nbparam->vdwtype))
862 ocl_copy_H2D_async(d_atdat->lj_comb, nbat->params().lj_comb.data(), 0, natoms * sizeof(cl_float2),
863 ls, bDoTime ? timers->atdat.fetchNextEvent() : nullptr);
865 else
867 ocl_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(), 0, natoms * sizeof(int),
868 ls, bDoTime ? timers->atdat.fetchNextEvent() : nullptr);
871 if (bDoTime)
873 timers->atdat.closeTimingRegion(ls);
876 /* kick off the tasks enqueued above to ensure concurrency with the search */
877 cl_error = clFlush(ls);
878 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
879 ("clFlush failed: " + ocl_get_error_string(cl_error)).c_str());
882 /*! \brief Releases an OpenCL kernel pointer */
883 static void free_kernel(cl_kernel* kernel_ptr)
885 cl_int gmx_unused cl_error;
887 GMX_ASSERT(kernel_ptr, "Need a valid kernel pointer");
889 if (*kernel_ptr)
891 cl_error = clReleaseKernel(*kernel_ptr);
892 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
893 ("clReleaseKernel failed: " + ocl_get_error_string(cl_error)).c_str());
895 *kernel_ptr = nullptr;
899 /*! \brief Releases a list of OpenCL kernel pointers */
900 static void free_kernels(cl_kernel* kernels, int count)
902 int i;
904 for (i = 0; i < count; i++)
906 free_kernel(kernels + i);
910 /*! \brief Free the OpenCL runtime data (context and program).
912 * The function releases the OpenCL context and program assuciated with the
913 * device that the calling PP rank is running on.
915 * \param runData [in] porinter to the structure with runtime data.
917 static void free_gpu_device_runtime_data(gmx_device_runtime_data_t* runData)
919 if (runData == nullptr)
921 return;
924 if (runData->program)
926 cl_int cl_error = clReleaseProgram(runData->program);
927 GMX_RELEASE_ASSERT(cl_error == CL_SUCCESS,
928 ("clReleaseProgram failed: " + ocl_get_error_string(cl_error)).c_str());
929 runData->program = nullptr;
933 //! This function is documented in the header file
934 void gpu_free(NbnxmGpu* nb)
936 if (nb == nullptr)
938 return;
941 /* Free kernels */
942 int kernel_count = sizeof(nb->kernel_ener_noprune_ptr) / sizeof(nb->kernel_ener_noprune_ptr[0][0]);
943 free_kernels(nb->kernel_ener_noprune_ptr[0], kernel_count);
945 kernel_count = sizeof(nb->kernel_ener_prune_ptr) / sizeof(nb->kernel_ener_prune_ptr[0][0]);
946 free_kernels(nb->kernel_ener_prune_ptr[0], kernel_count);
948 kernel_count = sizeof(nb->kernel_noener_noprune_ptr) / sizeof(nb->kernel_noener_noprune_ptr[0][0]);
949 free_kernels(nb->kernel_noener_noprune_ptr[0], kernel_count);
951 kernel_count = sizeof(nb->kernel_noener_prune_ptr) / sizeof(nb->kernel_noener_prune_ptr[0][0]);
952 free_kernels(nb->kernel_noener_prune_ptr[0], kernel_count);
954 free_kernel(&(nb->kernel_zero_e_fshift));
956 /* Free atdat */
957 freeDeviceBuffer(&(nb->atdat->xq));
958 freeDeviceBuffer(&(nb->atdat->f));
959 freeDeviceBuffer(&(nb->atdat->e_lj));
960 freeDeviceBuffer(&(nb->atdat->e_el));
961 freeDeviceBuffer(&(nb->atdat->fshift));
962 freeDeviceBuffer(&(nb->atdat->lj_comb));
963 freeDeviceBuffer(&(nb->atdat->atom_types));
964 freeDeviceBuffer(&(nb->atdat->shift_vec));
965 sfree(nb->atdat);
967 /* Free nbparam */
968 freeDeviceBuffer(&(nb->nbparam->nbfp_climg2d));
969 freeDeviceBuffer(&(nb->nbparam->nbfp_comb_climg2d));
970 freeDeviceBuffer(&(nb->nbparam->coulomb_tab_climg2d));
971 sfree(nb->nbparam);
973 /* Free plist */
974 auto* plist = nb->plist[InteractionLocality::Local];
975 freeDeviceBuffer(&plist->sci);
976 freeDeviceBuffer(&plist->cj4);
977 freeDeviceBuffer(&plist->imask);
978 freeDeviceBuffer(&plist->excl);
979 sfree(plist);
980 if (nb->bUseTwoStreams)
982 auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
983 freeDeviceBuffer(&plist_nl->sci);
984 freeDeviceBuffer(&plist_nl->cj4);
985 freeDeviceBuffer(&plist_nl->imask);
986 freeDeviceBuffer(&plist_nl->excl);
987 sfree(plist_nl);
990 /* Free nbst */
991 pfree(nb->nbst.e_lj);
992 nb->nbst.e_lj = nullptr;
994 pfree(nb->nbst.e_el);
995 nb->nbst.e_el = nullptr;
997 pfree(nb->nbst.fshift);
998 nb->nbst.fshift = nullptr;
1000 /* Free command queues */
1001 clReleaseCommandQueue(nb->stream[InteractionLocality::Local]);
1002 nb->stream[InteractionLocality::Local] = nullptr;
1003 if (nb->bUseTwoStreams)
1005 clReleaseCommandQueue(nb->stream[InteractionLocality::NonLocal]);
1006 nb->stream[InteractionLocality::NonLocal] = nullptr;
1008 /* Free other events */
1009 if (nb->nonlocal_done)
1011 clReleaseEvent(nb->nonlocal_done);
1012 nb->nonlocal_done = nullptr;
1014 if (nb->misc_ops_and_local_H2D_done)
1016 clReleaseEvent(nb->misc_ops_and_local_H2D_done);
1017 nb->misc_ops_and_local_H2D_done = nullptr;
1020 free_gpu_device_runtime_data(nb->dev_rundata);
1021 delete nb->dev_rundata;
1023 /* Free timers and timings */
1024 delete nb->timers;
1025 sfree(nb->timings);
1026 delete nb;
1028 if (debug)
1030 fprintf(debug, "Cleaned up OpenCL data structures.\n");
1034 //! This function is documented in the header file
1035 gmx_wallclock_gpu_nbnxn_t* gpu_get_timings(NbnxmGpu* nb)
1037 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
1040 //! This function is documented in the header file
1041 void gpu_reset_timings(nonbonded_verlet_t* nbv)
1043 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
1045 init_timings(nbv->gpu_nbv->timings);
1049 //! This function is documented in the header file
1050 int gpu_min_ci_balanced(NbnxmGpu* nb)
1052 return nb != nullptr ? gpu_min_ci_balanced_factor * nb->deviceInfo->compute_units : 0;
1055 //! This function is documented in the header file
1056 gmx_bool gpu_is_kernel_ewald_analytical(const NbnxmGpu* nb)
1058 return ((nb->nbparam->eeltype == eelOclEWALD_ANA) || (nb->nbparam->eeltype == eelOclEWALD_ANA_TWIN));
1061 } // namespace Nbnxm