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.
37 * \brief Define CUDA implementation of nbnxn_gpu_data_mgmt.h
39 * \author Szilard Pall <pall.szilard@gmail.com>
48 // TODO We would like to move this down, but the way NbnxmGpu
49 // is currently declared means this has to be before gpu_types.h
50 #include "nbnxm_cuda_types.h"
52 // TODO Remove this comment when the above order issue is resolved
53 #include "gromacs/gpu_utils/cudautils.cuh"
54 #include "gromacs/gpu_utils/gpu_utils.h"
55 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
56 #include "gromacs/gpu_utils/pmalloc_cuda.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/gridset.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/basedefinitions.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/real.h"
74 #include "gromacs/utility/smalloc.h"
76 #include "nbnxm_cuda.h"
81 /* This is a heuristically determined parameter for the Kepler
82 * and Maxwell architectures for the minimum size of ci lists by multiplying
83 * this constant with the # of multiprocessors on the current device.
84 * Since the maximum number of blocks per multiprocessor is 16, the ideal
85 * count for small systems is 32 or 48 blocks per multiprocessor. Because
86 * there is a bit of fluctuations in the generated block counts, we use
87 * a target of 44 instead of the ideal value of 48.
89 static unsigned int gpu_min_ci_balanced_factor = 44;
92 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb);
95 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t* nbparam);
97 /*! \brief Return whether combination rules are used.
99 * \param[in] pointer to nonbonded paramter struct
100 * \return true if combination rules are used in this run, false otherwise
102 static inline bool useLjCombRule(const cu_nbparam_t* nbparam)
104 return (nbparam->vdwtype == evdwCuCUTCOMBGEOM || nbparam->vdwtype == evdwCuCUTCOMBLB);
107 /*! \brief Initialized the Ewald Coulomb correction GPU table.
109 Tabulates the Ewald Coulomb force and initializes the size/scale
110 and the table GPU array. If called with an already allocated table,
111 it just re-uploads the table.
113 static void init_ewald_coulomb_force_table(const EwaldCorrectionTables& tables, cu_nbparam_t* nbp)
115 if (nbp->coulomb_tab != nullptr)
117 nbnxn_cuda_free_nbparam_table(nbp);
120 nbp->coulomb_tab_scale = tables.scale;
121 initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj, tables.tableF.data(),
122 tables.tableF.size());
126 /*! Initializes the atomdata structure first time, it only gets filled at
128 static void init_atomdata_first(cu_atomdata_t* ad, int ntypes)
133 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS * sizeof(*ad->shift_vec));
134 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
135 ad->bShiftVecUploaded = false;
137 stat = cudaMalloc((void**)&ad->fshift, SHIFTS * sizeof(*ad->fshift));
138 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
140 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
141 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
142 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
143 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
145 /* initialize to nullptr poiters to data that is not allocated here and will
146 need reallocation in nbnxn_cuda_init_atomdata */
150 /* size -1 indicates that the respective array hasn't been initialized yet */
155 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
156 earlier GPUs, single or twin cut-off. */
157 static int pick_ewald_kernel_type(const interaction_const_t& ic)
159 bool bTwinCut = (ic.rcoulomb != ic.rvdw);
160 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
163 /* Benchmarking/development environment variables to force the use of
164 analytical or tabulated Ewald kernel. */
165 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != nullptr);
166 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != nullptr);
168 if (bForceAnalyticalEwald && bForceTabulatedEwald)
171 "Both analytical and tabulated Ewald CUDA non-bonded kernels "
172 "requested through environment variables.");
175 /* By default use analytical Ewald. */
176 bUseAnalyticalEwald = true;
177 if (bForceAnalyticalEwald)
181 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
184 else if (bForceTabulatedEwald)
186 bUseAnalyticalEwald = false;
190 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
194 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
195 forces it (use it for debugging/benchmarking only). */
196 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == nullptr))
198 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
202 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
208 /*! Copies all parameters related to the cut-off from ic to nbp */
209 static void set_cutoff_parameters(cu_nbparam_t* nbp, const interaction_const_t* ic, const PairlistParams& listParams)
211 nbp->ewald_beta = ic->ewaldcoeff_q;
212 nbp->sh_ewald = ic->sh_ewald;
213 nbp->epsfac = ic->epsfac;
214 nbp->two_k_rf = 2.0 * ic->k_rf;
215 nbp->c_rf = ic->c_rf;
216 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
217 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
218 nbp->rlistOuter_sq = listParams.rlistOuter * listParams.rlistOuter;
219 nbp->rlistInner_sq = listParams.rlistInner * listParams.rlistInner;
220 nbp->useDynamicPruning = listParams.useDynamicPruning;
222 nbp->sh_lj_ewald = ic->sh_lj_ewald;
223 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
225 nbp->rvdw_switch = ic->rvdw_switch;
226 nbp->dispersion_shift = ic->dispersion_shift;
227 nbp->repulsion_shift = ic->repulsion_shift;
228 nbp->vdw_switch = ic->vdw_switch;
231 /*! Initializes the nonbonded parameter data structure. */
232 static void init_nbparam(cu_nbparam_t* nbp,
233 const interaction_const_t* ic,
234 const PairlistParams& listParams,
235 const nbnxn_atomdata_t::Params& nbatParams)
239 ntypes = nbatParams.numTypes;
241 set_cutoff_parameters(nbp, ic, listParams);
243 /* The kernel code supports LJ combination rules (geometric and LB) for
244 * all kernel types, but we only generate useful combination rule kernels.
245 * We currently only use LJ combination rule (geometric and LB) kernels
246 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
247 * with PME and 20% with RF, the other kernels speed up about half as much.
248 * For LJ force-switch the geometric rule would give 7% speed-up, but this
249 * combination is rarely used. LJ force-switch with LB rule is more common,
250 * but gives only 1% speed-up.
252 if (ic->vdwtype == evdwCUT)
254 switch (ic->vdw_modifier)
257 case eintmodPOTSHIFT:
258 switch (nbatParams.comb_rule)
260 case ljcrNONE: nbp->vdwtype = evdwCuCUT; break;
261 case ljcrGEOM: nbp->vdwtype = evdwCuCUTCOMBGEOM; break;
262 case ljcrLB: nbp->vdwtype = evdwCuCUTCOMBLB; break;
265 "The requested LJ combination rule is not implemented in the CUDA "
266 "GPU accelerated kernels!");
269 case eintmodFORCESWITCH: nbp->vdwtype = evdwCuFSWITCH; break;
270 case eintmodPOTSWITCH: nbp->vdwtype = evdwCuPSWITCH; break;
273 "The requested VdW interaction modifier is not implemented in the CUDA GPU "
274 "accelerated kernels!");
277 else if (ic->vdwtype == evdwPME)
279 if (ic->ljpme_comb_rule == ljcrGEOM)
281 assert(nbatParams.comb_rule == ljcrGEOM);
282 nbp->vdwtype = evdwCuEWALDGEOM;
286 assert(nbatParams.comb_rule == ljcrLB);
287 nbp->vdwtype = evdwCuEWALDLB;
293 "The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
296 if (ic->eeltype == eelCUT)
298 nbp->eeltype = eelCuCUT;
300 else if (EEL_RF(ic->eeltype))
302 nbp->eeltype = eelCuRF;
304 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
306 nbp->eeltype = pick_ewald_kernel_type(*ic);
310 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
312 "The requested electrostatics type is not implemented in the CUDA GPU accelerated "
316 /* generate table for PME */
317 nbp->coulomb_tab = nullptr;
318 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
320 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
321 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp);
324 /* set up LJ parameter lookup table */
325 if (!useLjCombRule(nbp))
327 initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj, nbatParams.nbfp.data(), 2 * ntypes * ntypes);
330 /* set up LJ-PME parameter lookup table */
331 if (ic->vdwtype == evdwPME)
333 initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj, nbatParams.nbfp_comb.data(), 2 * ntypes);
337 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
338 * electrostatic type switching to twin cut-off (or back) if needed. */
339 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t* nbv, const interaction_const_t* ic)
341 if (!nbv || !nbv->useGpu())
345 cu_nbparam_t* nbp = nbv->gpu_nbv->nbparam;
347 set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
349 nbp->eeltype = pick_ewald_kernel_type(*ic);
351 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
352 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp);
355 /*! Initializes the pair list data structure. */
356 static void init_plist(cu_plist_t* pl)
358 /* initialize to nullptr pointers to data that is not allocated here and will
359 need reallocation in nbnxn_gpu_init_pairlist */
365 /* size -1 indicates that the respective array hasn't been initialized yet */
372 pl->imask_nalloc = -1;
374 pl->excl_nalloc = -1;
375 pl->haveFreshList = false;
378 /*! Initializes the timings data structure. */
379 static void init_timings(gmx_wallclock_gpu_nbnxn_t* t)
388 for (i = 0; i < 2; i++)
390 for (j = 0; j < 2; j++)
392 t->ktime[i][j].t = 0.0;
393 t->ktime[i][j].c = 0;
397 t->pruneTime.t = 0.0;
398 t->dynamicPruneTime.c = 0;
399 t->dynamicPruneTime.t = 0.0;
402 /*! Initializes simulation constant data. */
403 static void cuda_init_const(NbnxmGpu* nb,
404 const interaction_const_t* ic,
405 const PairlistParams& listParams,
406 const nbnxn_atomdata_t::Params& nbatParams)
408 init_atomdata_first(nb->atdat, nbatParams.numTypes);
409 init_nbparam(nb->nbparam, ic, listParams, nbatParams);
411 /* clear energy and shift force outputs */
412 nbnxn_cuda_clear_e_fshift(nb);
415 NbnxmGpu* gpu_init(const gmx_device_info_t* deviceInfo,
416 const interaction_const_t* ic,
417 const PairlistParams& listParams,
418 const nbnxn_atomdata_t* nbat,
420 bool bLocalAndNonlocal)
424 auto nb = new NbnxmGpu;
426 snew(nb->nbparam, 1);
427 snew(nb->plist[InteractionLocality::Local], 1);
428 if (bLocalAndNonlocal)
430 snew(nb->plist[InteractionLocality::NonLocal], 1);
433 nb->bUseTwoStreams = bLocalAndNonlocal;
435 nb->timers = new cu_timers_t();
436 snew(nb->timings, 1);
439 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
440 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
441 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
443 init_plist(nb->plist[InteractionLocality::Local]);
445 /* set device info, just point it to the right GPU among the detected ones */
446 nb->dev_info = deviceInfo;
448 /* local/non-local GPU streams */
449 stat = cudaStreamCreate(&nb->stream[InteractionLocality::Local]);
450 CU_RET_ERR(stat, "cudaStreamCreate on stream[InterationLocality::Local] failed");
451 if (nb->bUseTwoStreams)
453 init_plist(nb->plist[InteractionLocality::NonLocal]);
455 /* Note that the device we're running on does not have to support
456 * priorities, because we are querying the priority range which in this
457 * case will be a single value.
459 int highest_priority;
460 stat = cudaDeviceGetStreamPriorityRange(nullptr, &highest_priority);
461 CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
463 stat = cudaStreamCreateWithPriority(&nb->stream[InteractionLocality::NonLocal],
464 cudaStreamDefault, highest_priority);
466 "cudaStreamCreateWithPriority on stream[InteractionLocality::NonLocal] failed");
469 /* init events for sychronization (timing disabled for performance reasons!) */
470 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
471 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
472 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
473 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
475 nb->xNonLocalCopyD2HDone = new GpuEventSynchronizer();
477 /* WARNING: CUDA timings are incorrect with multiple streams.
478 * This is the main reason why they are disabled by default.
480 // TODO: Consider turning on by default when we can detect nr of streams.
481 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
485 init_timings(nb->timings);
488 /* set the kernel type for the current GPU */
489 /* pick L1 cache configuration */
490 cuda_set_cacheconfig();
492 cuda_init_const(nb, ic, listParams, nbat->params());
494 nb->atomIndicesSize = 0;
495 nb->atomIndicesSize_alloc = 0;
497 nb->ncxy_na_alloc = 0;
499 nb->ncxy_ind_alloc = 0;
505 fprintf(debug, "Initialized CUDA data structures.\n");
511 void gpu_init_pairlist(NbnxmGpu* nb, const NbnxnPairlistGpu* h_plist, const InteractionLocality iloc)
514 bool bDoTime = (nb->bDoTime && !h_plist->sci.empty());
515 cudaStream_t stream = nb->stream[iloc];
516 cu_plist_t* d_plist = nb->plist[iloc];
518 if (d_plist->na_c < 0)
520 d_plist->na_c = h_plist->na_ci;
524 if (d_plist->na_c != h_plist->na_ci)
526 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
527 d_plist->na_c, h_plist->na_ci);
532 gpu_timers_t::Interaction& iTimers = nb->timers->interaction[iloc];
536 iTimers.pl_h2d.openTimingRegion(stream);
537 iTimers.didPairlistH2D = true;
540 DeviceContext context = nullptr;
542 reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(), &d_plist->nsci, &d_plist->sci_nalloc, context);
543 copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(), stream,
544 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
546 reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(), &d_plist->ncj4, &d_plist->cj4_nalloc, context);
547 copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(), stream,
548 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
550 reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size() * c_nbnxnGpuClusterpairSplit,
551 &d_plist->nimask, &d_plist->imask_nalloc, context);
553 reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(), &d_plist->nexcl,
554 &d_plist->excl_nalloc, context);
555 copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(), stream,
556 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
560 iTimers.pl_h2d.closeTimingRegion(stream);
563 /* the next use of thist list we be the first one, so we need to prune */
564 d_plist->haveFreshList = true;
567 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
569 cu_atomdata_t* adat = nb->atdat;
570 cudaStream_t ls = nb->stream[InteractionLocality::Local];
572 /* only if we have a dynamic box */
573 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
575 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(), SHIFTS * sizeof(*adat->shift_vec), ls);
576 adat->bShiftVecUploaded = true;
580 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
581 static void nbnxn_cuda_clear_f(NbnxmGpu* nb, int natoms_clear)
584 cu_atomdata_t* adat = nb->atdat;
585 cudaStream_t ls = nb->stream[InteractionLocality::Local];
587 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
588 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
591 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
592 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb)
595 cu_atomdata_t* adat = nb->atdat;
596 cudaStream_t ls = nb->stream[InteractionLocality::Local];
598 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
599 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
600 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
601 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
602 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
603 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
606 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
608 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
609 /* clear shift force array and energies if the outputs were
610 used in the current step */
613 nbnxn_cuda_clear_e_fshift(nb);
617 void gpu_init_atomdata(NbnxmGpu* nb, const nbnxn_atomdata_t* nbat)
622 bool bDoTime = nb->bDoTime;
623 cu_timers_t* timers = nb->timers;
624 cu_atomdata_t* d_atdat = nb->atdat;
625 cudaStream_t ls = nb->stream[InteractionLocality::Local];
627 natoms = nbat->numAtoms();
632 /* time async copy */
633 timers->atdat.openTimingRegion(ls);
636 /* need to reallocate if we have to copy more atoms than the amount of space
637 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
638 if (natoms > d_atdat->nalloc)
640 nalloc = over_alloc_small(natoms);
642 /* free up first if the arrays have already been initialized */
643 if (d_atdat->nalloc != -1)
645 freeDeviceBuffer(&d_atdat->f);
646 freeDeviceBuffer(&d_atdat->xq);
647 freeDeviceBuffer(&d_atdat->atom_types);
648 freeDeviceBuffer(&d_atdat->lj_comb);
651 stat = cudaMalloc((void**)&d_atdat->f, nalloc * sizeof(*d_atdat->f));
652 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
653 stat = cudaMalloc((void**)&d_atdat->xq, nalloc * sizeof(*d_atdat->xq));
654 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
655 if (useLjCombRule(nb->nbparam))
657 stat = cudaMalloc((void**)&d_atdat->lj_comb, nalloc * sizeof(*d_atdat->lj_comb));
658 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
662 stat = cudaMalloc((void**)&d_atdat->atom_types, nalloc * sizeof(*d_atdat->atom_types));
663 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
666 d_atdat->nalloc = nalloc;
670 d_atdat->natoms = natoms;
671 d_atdat->natoms_local = nbat->natoms_local;
673 /* need to clear GPU f output if realloc happened */
676 nbnxn_cuda_clear_f(nb, nalloc);
679 if (useLjCombRule(nb->nbparam))
681 cu_copy_H2D_async(d_atdat->lj_comb, nbat->params().lj_comb.data(),
682 natoms * sizeof(*d_atdat->lj_comb), ls);
686 cu_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(),
687 natoms * sizeof(*d_atdat->atom_types), ls);
692 timers->atdat.closeTimingRegion(ls);
696 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t* nbparam)
698 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
700 destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
704 void gpu_free(NbnxmGpu* nb)
707 cu_atomdata_t* atdat;
708 cu_nbparam_t* nbparam;
716 nbparam = nb->nbparam;
718 nbnxn_cuda_free_nbparam_table(nbparam);
720 stat = cudaEventDestroy(nb->nonlocal_done);
721 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
722 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
723 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
728 /* The non-local counters/stream (second in the array) are needed only with DD. */
729 for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
731 stat = cudaStreamDestroy(nb->stream[i]);
732 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
736 if (!useLjCombRule(nb->nbparam))
738 destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj);
741 if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
743 destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
746 stat = cudaFree(atdat->shift_vec);
747 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
748 stat = cudaFree(atdat->fshift);
749 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
751 stat = cudaFree(atdat->e_lj);
752 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
753 stat = cudaFree(atdat->e_el);
754 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
756 freeDeviceBuffer(&atdat->f);
757 freeDeviceBuffer(&atdat->xq);
758 freeDeviceBuffer(&atdat->atom_types);
759 freeDeviceBuffer(&atdat->lj_comb);
762 auto* plist = nb->plist[InteractionLocality::Local];
763 freeDeviceBuffer(&plist->sci);
764 freeDeviceBuffer(&plist->cj4);
765 freeDeviceBuffer(&plist->imask);
766 freeDeviceBuffer(&plist->excl);
768 if (nb->bUseTwoStreams)
770 auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
771 freeDeviceBuffer(&plist_nl->sci);
772 freeDeviceBuffer(&plist_nl->cj4);
773 freeDeviceBuffer(&plist_nl->imask);
774 freeDeviceBuffer(&plist_nl->excl);
779 pfree(nb->nbst.e_lj);
780 nb->nbst.e_lj = nullptr;
782 pfree(nb->nbst.e_el);
783 nb->nbst.e_el = nullptr;
785 pfree(nb->nbst.fshift);
786 nb->nbst.fshift = nullptr;
795 fprintf(debug, "Cleaned up CUDA data structures.\n");
799 //! This function is documented in the header file
800 gmx_wallclock_gpu_nbnxn_t* gpu_get_timings(NbnxmGpu* nb)
802 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
805 void gpu_reset_timings(nonbonded_verlet_t* nbv)
807 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
809 init_timings(nbv->gpu_nbv->timings);
813 int gpu_min_ci_balanced(NbnxmGpu* nb)
815 return nb != nullptr ? gpu_min_ci_balanced_factor * nb->dev_info->prop.multiProcessorCount : 0;
818 gmx_bool gpu_is_kernel_ewald_analytical(const NbnxmGpu* nb)
820 return ((nb->nbparam->eeltype == eelCuEWALD_ANA) || (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
823 void* gpu_get_command_stream(NbnxmGpu* nb, const InteractionLocality iloc)
827 return static_cast<void*>(&nb->stream[iloc]);
830 void* gpu_get_xq(NbnxmGpu* nb)
834 return static_cast<void*>(nb->atdat->xq);
837 void* gpu_get_f(NbnxmGpu* nb)
841 return static_cast<void*>(nb->atdat->f);
844 rvec* gpu_get_fshift(NbnxmGpu* nb)
848 return reinterpret_cast<rvec*>(nb->atdat->fshift);
851 /* Initialization for X buffer operations on GPU. */
852 /* TODO Remove explicit pinning from host arrays from here and manage in a more natural way*/
853 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet& gridSet, NbnxmGpu* gpu_nbv)
855 cudaStream_t stream = gpu_nbv->stream[InteractionLocality::Local];
856 bool bDoTime = gpu_nbv->bDoTime;
857 const int maxNumColumns = gridSet.numColumnsMax();
859 reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns * gridSet.grids().size(),
860 &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, nullptr);
861 reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns * gridSet.grids().size(),
862 &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, nullptr);
864 for (unsigned int g = 0; g < gridSet.grids().size(); g++)
867 const Nbnxm::Grid& grid = gridSet.grids()[g];
869 const int numColumns = grid.numColumns();
870 const int* atomIndices = gridSet.atomIndices().data();
871 const int atomIndicesSize = gridSet.atomIndices().size();
872 const int* cxy_na = grid.cxy_na().data();
873 const int* cxy_ind = grid.cxy_ind().data();
875 reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize,
876 &gpu_nbv->atomIndicesSize_alloc, nullptr);
878 if (atomIndicesSize > 0)
883 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
886 copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, stream,
887 GpuApiCallBehavior::Async, nullptr);
891 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
899 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
902 int* destPtr = &gpu_nbv->cxy_na[maxNumColumns * g];
903 copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, stream, GpuApiCallBehavior::Async, nullptr);
907 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
912 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
915 destPtr = &gpu_nbv->cxy_ind[maxNumColumns * g];
916 copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, stream, GpuApiCallBehavior::Async, nullptr);
920 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
925 // The above data is transferred on the local stream but is a
926 // dependency of the nonlocal stream (specifically the nonlocal X
927 // buf ops kernel). We therefore set a dependency to ensure
928 // that the nonlocal stream waits on the local stream here.
929 // This call records an event in the local stream:
930 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
931 // ...and this call instructs the nonlocal stream to wait on that event:
932 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
937 /* Initialization for F buffer operations on GPU. */
938 void nbnxn_gpu_init_add_nbat_f_to_f(const int* cell,
941 GpuEventSynchronizer* const localReductionDone)
944 cudaStream_t stream = gpu_nbv->stream[InteractionLocality::Local];
946 GMX_ASSERT(localReductionDone, "localReductionDone should be a valid pointer");
947 gpu_nbv->localFReductionDone = localReductionDone;
949 if (natoms_total > 0)
951 reallocateDeviceBuffer(&gpu_nbv->cell, natoms_total, &gpu_nbv->ncell, &gpu_nbv->ncell_alloc, nullptr);
952 copyToDeviceBuffer(&gpu_nbv->cell, cell, 0, natoms_total, stream, GpuApiCallBehavior::Async, nullptr);