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/device_context.h"
55 #include "gromacs/gpu_utils/device_stream_manager.h"
56 #include "gromacs/gpu_utils/gpu_utils.h"
57 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
58 #include "gromacs/gpu_utils/pmalloc_cuda.h"
59 #include "gromacs/hardware/device_information.h"
60 #include "gromacs/hardware/device_management.h"
61 #include "gromacs/math/vectypes.h"
62 #include "gromacs/mdlib/force_flags.h"
63 #include "gromacs/mdtypes/interaction_const.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/nbnxm/atomdata.h"
66 #include "gromacs/nbnxm/gpu_data_mgmt.h"
67 #include "gromacs/nbnxm/gridset.h"
68 #include "gromacs/nbnxm/nbnxm.h"
69 #include "gromacs/nbnxm/nbnxm_gpu.h"
70 #include "gromacs/nbnxm/nbnxm_gpu_data_mgmt.h"
71 #include "gromacs/nbnxm/pairlistsets.h"
72 #include "gromacs/pbcutil/ishift.h"
73 #include "gromacs/timing/gpu_timing.h"
74 #include "gromacs/utility/basedefinitions.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/real.h"
78 #include "gromacs/utility/smalloc.h"
80 #include "nbnxm_cuda.h"
85 /* This is a heuristically determined parameter for the Kepler
86 * and Maxwell architectures for the minimum size of ci lists by multiplying
87 * this constant with the # of multiprocessors on the current device.
88 * Since the maximum number of blocks per multiprocessor is 16, the ideal
89 * count for small systems is 32 or 48 blocks per multiprocessor. Because
90 * there is a bit of fluctuations in the generated block counts, we use
91 * a target of 44 instead of the ideal value of 48.
93 static unsigned int gpu_min_ci_balanced_factor = 44;
96 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb);
98 /*! Initializes the atomdata structure first time, it only gets filled at
100 static void init_atomdata_first(cu_atomdata_t* ad, int ntypes, const DeviceContext& deviceContext)
103 allocateDeviceBuffer(&ad->shift_vec, SHIFTS, deviceContext);
104 ad->bShiftVecUploaded = false;
106 allocateDeviceBuffer(&ad->fshift, SHIFTS, deviceContext);
107 allocateDeviceBuffer(&ad->e_lj, 1, deviceContext);
108 allocateDeviceBuffer(&ad->e_el, 1, deviceContext);
110 /* initialize to nullptr poiters to data that is not allocated here and will
111 need reallocation in nbnxn_cuda_init_atomdata */
115 /* size -1 indicates that the respective array hasn't been initialized yet */
120 /*! Initializes the nonbonded parameter data structure. */
121 static void init_nbparam(NBParamGpu* nbp,
122 const interaction_const_t* ic,
123 const PairlistParams& listParams,
124 const nbnxn_atomdata_t::Params& nbatParams,
125 const DeviceContext& deviceContext)
127 const int ntypes = nbatParams.numTypes;
129 set_cutoff_parameters(nbp, ic, listParams);
131 /* The kernel code supports LJ combination rules (geometric and LB) for
132 * all kernel types, but we only generate useful combination rule kernels.
133 * We currently only use LJ combination rule (geometric and LB) kernels
134 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
135 * with PME and 20% with RF, the other kernels speed up about half as much.
136 * For LJ force-switch the geometric rule would give 7% speed-up, but this
137 * combination is rarely used. LJ force-switch with LB rule is more common,
138 * but gives only 1% speed-up.
140 nbp->vdwType = nbnxmGpuPickVdwKernelType(ic, nbatParams.comb_rule);
141 nbp->elecType = nbnxmGpuPickElectrostaticsKernelType(ic);
143 /* generate table for PME */
144 nbp->coulomb_tab = nullptr;
145 if (nbp->elecType == ElecType::EwaldTab || nbp->elecType == ElecType::EwaldTabTwin)
147 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
148 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, deviceContext);
151 /* set up LJ parameter lookup table */
152 if (!useLjCombRule(nbp->vdwType))
154 initParamLookupTable(&nbp->nbfp, &nbp->nbfp_texobj, nbatParams.nbfp.data(),
155 2 * ntypes * ntypes, deviceContext);
158 /* set up LJ-PME parameter lookup table */
159 if (ic->vdwtype == evdwPME)
161 initParamLookupTable(&nbp->nbfp_comb, &nbp->nbfp_comb_texobj, nbatParams.nbfp_comb.data(),
162 2 * ntypes, deviceContext);
166 /*! Initializes simulation constant data. */
167 static void cuda_init_const(NbnxmGpu* nb,
168 const interaction_const_t* ic,
169 const PairlistParams& listParams,
170 const nbnxn_atomdata_t::Params& nbatParams)
172 init_atomdata_first(nb->atdat, nbatParams.numTypes, *nb->deviceContext_);
173 init_nbparam(nb->nbparam, ic, listParams, nbatParams, *nb->deviceContext_);
175 /* clear energy and shift force outputs */
176 nbnxn_cuda_clear_e_fshift(nb);
179 NbnxmGpu* gpu_init(const gmx::DeviceStreamManager& deviceStreamManager,
180 const interaction_const_t* ic,
181 const PairlistParams& listParams,
182 const nbnxn_atomdata_t* nbat,
183 bool bLocalAndNonlocal)
187 auto nb = new NbnxmGpu();
188 nb->deviceContext_ = &deviceStreamManager.context();
190 snew(nb->nbparam, 1);
191 snew(nb->plist[InteractionLocality::Local], 1);
192 if (bLocalAndNonlocal)
194 snew(nb->plist[InteractionLocality::NonLocal], 1);
197 nb->bUseTwoStreams = bLocalAndNonlocal;
199 nb->timers = new cu_timers_t();
200 snew(nb->timings, 1);
203 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
204 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
205 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
207 init_plist(nb->plist[InteractionLocality::Local]);
209 /* local/non-local GPU streams */
210 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedLocal),
211 "Local non-bonded stream should be initialized to use GPU for non-bonded.");
212 nb->deviceStreams[InteractionLocality::Local] =
213 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedLocal);
214 if (nb->bUseTwoStreams)
216 init_plist(nb->plist[InteractionLocality::NonLocal]);
218 /* Note that the device we're running on does not have to support
219 * priorities, because we are querying the priority range which in this
220 * case will be a single value.
222 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedNonLocal),
223 "Non-local non-bonded stream should be initialized to use GPU for "
224 "non-bonded with domain decomposition.");
225 nb->deviceStreams[InteractionLocality::NonLocal] =
226 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedNonLocal);
230 /* init events for sychronization (timing disabled for performance reasons!) */
231 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
232 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
233 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
234 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
236 nb->xNonLocalCopyD2HDone = new GpuEventSynchronizer();
238 /* WARNING: CUDA timings are incorrect with multiple streams.
239 * This is the main reason why they are disabled by default.
241 // TODO: Consider turning on by default when we can detect nr of streams.
242 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
246 init_timings(nb->timings);
249 /* set the kernel type for the current GPU */
250 /* pick L1 cache configuration */
251 cuda_set_cacheconfig();
253 cuda_init_const(nb, ic, listParams, nbat->params());
255 nb->atomIndicesSize = 0;
256 nb->atomIndicesSize_alloc = 0;
258 nb->ncxy_na_alloc = 0;
260 nb->ncxy_ind_alloc = 0;
264 fprintf(debug, "Initialized CUDA data structures.\n");
270 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
272 cu_atomdata_t* adat = nb->atdat;
273 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
275 /* only if we have a dynamic box */
276 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
278 static_assert(sizeof(adat->shift_vec[0]) == sizeof(nbatom->shift_vec[0]),
279 "Sizes of host- and device-side shift vectors should be the same.");
280 copyToDeviceBuffer(&adat->shift_vec, reinterpret_cast<const float3*>(nbatom->shift_vec.data()),
281 0, SHIFTS, localStream, GpuApiCallBehavior::Async, nullptr);
282 adat->bShiftVecUploaded = true;
286 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
287 static void nbnxn_cuda_clear_f(NbnxmGpu* nb, int natoms_clear)
289 cu_atomdata_t* adat = nb->atdat;
290 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
291 clearDeviceBufferAsync(&adat->f, 0, natoms_clear, localStream);
294 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
295 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb)
297 cu_atomdata_t* adat = nb->atdat;
298 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
300 clearDeviceBufferAsync(&adat->fshift, 0, SHIFTS, localStream);
301 clearDeviceBufferAsync(&adat->e_lj, 0, 1, localStream);
302 clearDeviceBufferAsync(&adat->e_el, 0, 1, localStream);
305 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
307 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
308 /* clear shift force array and energies if the outputs were
309 used in the current step */
312 nbnxn_cuda_clear_e_fshift(nb);
316 void gpu_init_atomdata(NbnxmGpu* nb, const nbnxn_atomdata_t* nbat)
320 bool bDoTime = nb->bDoTime;
321 cu_timers_t* timers = nb->timers;
322 cu_atomdata_t* d_atdat = nb->atdat;
323 const DeviceContext& deviceContext = *nb->deviceContext_;
324 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
326 natoms = nbat->numAtoms();
331 /* time async copy */
332 timers->atdat.openTimingRegion(localStream);
335 /* need to reallocate if we have to copy more atoms than the amount of space
336 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
337 if (natoms > d_atdat->nalloc)
339 nalloc = over_alloc_small(natoms);
341 /* free up first if the arrays have already been initialized */
342 if (d_atdat->nalloc != -1)
344 freeDeviceBuffer(&d_atdat->f);
345 freeDeviceBuffer(&d_atdat->xq);
346 freeDeviceBuffer(&d_atdat->atom_types);
347 freeDeviceBuffer(&d_atdat->lj_comb);
350 allocateDeviceBuffer(&d_atdat->f, nalloc, deviceContext);
351 allocateDeviceBuffer(&d_atdat->xq, nalloc, deviceContext);
352 if (useLjCombRule(nb->nbparam->vdwType))
354 allocateDeviceBuffer(&d_atdat->lj_comb, nalloc, deviceContext);
358 allocateDeviceBuffer(&d_atdat->atom_types, nalloc, deviceContext);
361 d_atdat->nalloc = nalloc;
365 d_atdat->natoms = natoms;
366 d_atdat->natoms_local = nbat->natoms_local;
368 /* need to clear GPU f output if realloc happened */
371 nbnxn_cuda_clear_f(nb, nalloc);
374 if (useLjCombRule(nb->nbparam->vdwType))
376 static_assert(sizeof(d_atdat->lj_comb[0]) == sizeof(float2),
377 "Size of the LJ parameters element should be equal to the size of float2.");
378 copyToDeviceBuffer(&d_atdat->lj_comb,
379 reinterpret_cast<const float2*>(nbat->params().lj_comb.data()), 0,
380 natoms, localStream, GpuApiCallBehavior::Async, nullptr);
384 static_assert(sizeof(d_atdat->atom_types[0]) == sizeof(nbat->params().type[0]),
385 "Sizes of host- and device-side atom types should be the same.");
386 copyToDeviceBuffer(&d_atdat->atom_types, nbat->params().type.data(), 0, natoms, localStream,
387 GpuApiCallBehavior::Async, nullptr);
392 timers->atdat.closeTimingRegion(localStream);
396 void gpu_free(NbnxmGpu* nb)
399 cu_atomdata_t* atdat;
408 nbparam = nb->nbparam;
410 if ((!nbparam->coulomb_tab)
411 && (nbparam->elecType == ElecType::EwaldTab || nbparam->elecType == ElecType::EwaldTabTwin))
413 destroyParamLookupTable(&nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
416 stat = cudaEventDestroy(nb->nonlocal_done);
417 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
418 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
419 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
423 if (!useLjCombRule(nb->nbparam->vdwType))
425 destroyParamLookupTable(&nbparam->nbfp, nbparam->nbfp_texobj);
428 if (nbparam->vdwType == VdwType::EwaldGeom || nbparam->vdwType == VdwType::EwaldLB)
430 destroyParamLookupTable(&nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
433 freeDeviceBuffer(&atdat->shift_vec);
434 freeDeviceBuffer(&atdat->fshift);
436 freeDeviceBuffer(&atdat->e_lj);
437 freeDeviceBuffer(&atdat->e_el);
439 freeDeviceBuffer(&atdat->f);
440 freeDeviceBuffer(&atdat->xq);
441 freeDeviceBuffer(&atdat->atom_types);
442 freeDeviceBuffer(&atdat->lj_comb);
445 auto* plist = nb->plist[InteractionLocality::Local];
446 freeDeviceBuffer(&plist->sci);
447 freeDeviceBuffer(&plist->cj4);
448 freeDeviceBuffer(&plist->imask);
449 freeDeviceBuffer(&plist->excl);
451 if (nb->bUseTwoStreams)
453 auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
454 freeDeviceBuffer(&plist_nl->sci);
455 freeDeviceBuffer(&plist_nl->cj4);
456 freeDeviceBuffer(&plist_nl->imask);
457 freeDeviceBuffer(&plist_nl->excl);
462 pfree(nb->nbst.e_lj);
463 nb->nbst.e_lj = nullptr;
465 pfree(nb->nbst.e_el);
466 nb->nbst.e_el = nullptr;
468 pfree(nb->nbst.fshift);
469 nb->nbst.fshift = nullptr;
478 fprintf(debug, "Cleaned up CUDA data structures.\n");
482 int gpu_min_ci_balanced(NbnxmGpu* nb)
484 return nb != nullptr ? gpu_min_ci_balanced_factor * nb->deviceContext_->deviceInfo().prop.multiProcessorCount
488 void* gpu_get_xq(NbnxmGpu* nb)
492 return static_cast<void*>(nb->atdat->xq);
495 DeviceBuffer<gmx::RVec> gpu_get_f(NbnxmGpu* nb)
499 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->f);
502 DeviceBuffer<gmx::RVec> gpu_get_fshift(NbnxmGpu* nb)
506 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->fshift);
509 /* Initialization for X buffer operations on GPU. */
510 /* TODO Remove explicit pinning from host arrays from here and manage in a more natural way*/
511 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet& gridSet, NbnxmGpu* gpu_nbv)
513 const DeviceStream& deviceStream = *gpu_nbv->deviceStreams[InteractionLocality::Local];
514 bool bDoTime = gpu_nbv->bDoTime;
515 const int maxNumColumns = gridSet.numColumnsMax();
517 reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns * gridSet.grids().size(),
518 &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, *gpu_nbv->deviceContext_);
519 reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns * gridSet.grids().size(),
520 &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, *gpu_nbv->deviceContext_);
522 for (unsigned int g = 0; g < gridSet.grids().size(); g++)
525 const Nbnxm::Grid& grid = gridSet.grids()[g];
527 const int numColumns = grid.numColumns();
528 const int* atomIndices = gridSet.atomIndices().data();
529 const int atomIndicesSize = gridSet.atomIndices().size();
530 const int* cxy_na = grid.cxy_na().data();
531 const int* cxy_ind = grid.cxy_ind().data();
533 reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize,
534 &gpu_nbv->atomIndicesSize_alloc, *gpu_nbv->deviceContext_);
536 if (atomIndicesSize > 0)
541 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
544 copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, deviceStream,
545 GpuApiCallBehavior::Async, nullptr);
549 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
557 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
560 int* destPtr = &gpu_nbv->cxy_na[maxNumColumns * g];
561 copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, deviceStream,
562 GpuApiCallBehavior::Async, nullptr);
566 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
571 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
574 destPtr = &gpu_nbv->cxy_ind[maxNumColumns * g];
575 copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, deviceStream,
576 GpuApiCallBehavior::Async, nullptr);
580 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
585 // The above data is transferred on the local stream but is a
586 // dependency of the nonlocal stream (specifically the nonlocal X
587 // buf ops kernel). We therefore set a dependency to ensure
588 // that the nonlocal stream waits on the local stream here.
589 // This call records an event in the local stream:
590 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
591 // ...and this call instructs the nonlocal stream to wait on that event:
592 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);