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,2021, 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.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(NBAtomData* ad, int ntypes, const DeviceContext& deviceContext)
102 ad->numTypes = ntypes;
103 allocateDeviceBuffer(&ad->shiftVec, SHIFTS, deviceContext);
104 ad->shiftVecUploaded = false;
106 allocateDeviceBuffer(&ad->fShift, SHIFTS, deviceContext);
107 allocateDeviceBuffer(&ad->eLJ, 1, deviceContext);
108 allocateDeviceBuffer(&ad->eElec, 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 */
117 ad->numAtomsAlloc = -1;
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.ljCombinationRule);
141 nbp->elecType = nbnxmGpuPickElectrostaticsKernelType(ic, deviceContext.deviceInfo());
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 static_assert(sizeof(decltype(nbp->nbfp)) == 2 * sizeof(decltype(*nbatParams.nbfp.data())),
155 "Mismatch in the size of host / device data types");
156 initParamLookupTable(&nbp->nbfp,
158 reinterpret_cast<const Float2*>(nbatParams.nbfp.data()),
163 /* set up LJ-PME parameter lookup table */
164 if (ic->vdwtype == VanDerWaalsType::Pme)
166 static_assert(sizeof(decltype(nbp->nbfp_comb))
167 == 2 * sizeof(decltype(*nbatParams.nbfp_comb.data())),
168 "Mismatch in the size of host / device data types");
169 initParamLookupTable(&nbp->nbfp_comb,
170 &nbp->nbfp_comb_texobj,
171 reinterpret_cast<const Float2*>(nbatParams.nbfp_comb.data()),
177 /*! Initializes simulation constant data. */
178 static void cuda_init_const(NbnxmGpu* nb,
179 const interaction_const_t* ic,
180 const PairlistParams& listParams,
181 const nbnxn_atomdata_t::Params& nbatParams)
183 init_atomdata_first(nb->atdat, nbatParams.numTypes, *nb->deviceContext_);
184 init_nbparam(nb->nbparam, ic, listParams, nbatParams, *nb->deviceContext_);
186 /* clear energy and shift force outputs */
187 nbnxn_cuda_clear_e_fshift(nb);
190 NbnxmGpu* gpu_init(const gmx::DeviceStreamManager& deviceStreamManager,
191 const interaction_const_t* ic,
192 const PairlistParams& listParams,
193 const nbnxn_atomdata_t* nbat,
194 bool bLocalAndNonlocal)
196 auto nb = new NbnxmGpu();
197 nb->deviceContext_ = &deviceStreamManager.context();
199 snew(nb->nbparam, 1);
200 snew(nb->plist[InteractionLocality::Local], 1);
201 if (bLocalAndNonlocal)
203 snew(nb->plist[InteractionLocality::NonLocal], 1);
206 nb->bUseTwoStreams = bLocalAndNonlocal;
208 nb->timers = new Nbnxm::GpuTimers();
209 snew(nb->timings, 1);
212 pmalloc((void**)&nb->nbst.eLJ, sizeof(*nb->nbst.eLJ));
213 pmalloc((void**)&nb->nbst.eElec, sizeof(*nb->nbst.eElec));
214 pmalloc((void**)&nb->nbst.fShift, SHIFTS * sizeof(*nb->nbst.fShift));
216 init_plist(nb->plist[InteractionLocality::Local]);
218 /* local/non-local GPU streams */
219 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedLocal),
220 "Local non-bonded stream should be initialized to use GPU for non-bonded.");
221 nb->deviceStreams[InteractionLocality::Local] =
222 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedLocal);
223 if (nb->bUseTwoStreams)
225 init_plist(nb->plist[InteractionLocality::NonLocal]);
227 /* Note that the device we're running on does not have to support
228 * priorities, because we are querying the priority range which in this
229 * case will be a single value.
231 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedNonLocal),
232 "Non-local non-bonded stream should be initialized to use GPU for "
233 "non-bonded with domain decomposition.");
234 nb->deviceStreams[InteractionLocality::NonLocal] =
235 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedNonLocal);
239 /* WARNING: CUDA timings are incorrect with multiple streams.
240 * This is the main reason why they are disabled by default.
242 // TODO: Consider turning on by default when we can detect nr of streams.
243 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
247 init_timings(nb->timings);
250 /* set the kernel type for the current GPU */
251 /* pick L1 cache configuration */
252 cuda_set_cacheconfig();
254 cuda_init_const(nb, ic, listParams, nbat->params());
256 nb->atomIndicesSize = 0;
257 nb->atomIndicesSize_alloc = 0;
259 nb->ncxy_na_alloc = 0;
261 nb->ncxy_ind_alloc = 0;
265 fprintf(debug, "Initialized CUDA data structures.\n");
271 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
273 NBAtomData* adat = nb->atdat;
274 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
276 /* only if we have a dynamic box */
277 if (nbatom->bDynamicBox || !adat->shiftVecUploaded)
279 static_assert(sizeof(adat->shiftVec[0]) == sizeof(nbatom->shift_vec[0]),
280 "Sizes of host- and device-side shift vectors should be the same.");
281 copyToDeviceBuffer(&adat->shiftVec,
282 reinterpret_cast<const Float3*>(nbatom->shift_vec.data()),
286 GpuApiCallBehavior::Async,
288 adat->shiftVecUploaded = true;
292 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
293 static void nbnxn_cuda_clear_f(NbnxmGpu* nb, int natoms_clear)
295 NBAtomData* adat = nb->atdat;
296 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
297 clearDeviceBufferAsync(&adat->f, 0, natoms_clear, localStream);
300 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
301 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb)
303 NBAtomData* adat = nb->atdat;
304 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
306 clearDeviceBufferAsync(&adat->fShift, 0, SHIFTS, localStream);
307 clearDeviceBufferAsync(&adat->eLJ, 0, 1, localStream);
308 clearDeviceBufferAsync(&adat->eElec, 0, 1, localStream);
311 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
313 nbnxn_cuda_clear_f(nb, nb->atdat->numAtoms);
314 /* clear shift force array and energies if the outputs were
315 used in the current step */
318 nbnxn_cuda_clear_e_fshift(nb);
322 void gpu_free(NbnxmGpu* nb)
329 NBAtomData* atdat = nb->atdat;
330 NBParamGpu* nbparam = nb->nbparam;
332 if ((!nbparam->coulomb_tab)
333 && (nbparam->elecType == ElecType::EwaldTab || nbparam->elecType == ElecType::EwaldTabTwin))
335 destroyParamLookupTable(&nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
340 if (!useLjCombRule(nb->nbparam->vdwType))
342 destroyParamLookupTable(&nbparam->nbfp, nbparam->nbfp_texobj);
345 if (nbparam->vdwType == VdwType::EwaldGeom || nbparam->vdwType == VdwType::EwaldLB)
347 destroyParamLookupTable(&nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
350 freeDeviceBuffer(&atdat->shiftVec);
351 freeDeviceBuffer(&atdat->fShift);
353 freeDeviceBuffer(&atdat->eLJ);
354 freeDeviceBuffer(&atdat->eElec);
356 freeDeviceBuffer(&atdat->f);
357 freeDeviceBuffer(&atdat->xq);
358 freeDeviceBuffer(&atdat->atomTypes);
359 freeDeviceBuffer(&atdat->ljComb);
362 auto* plist = nb->plist[InteractionLocality::Local];
363 freeDeviceBuffer(&plist->sci);
364 freeDeviceBuffer(&plist->cj4);
365 freeDeviceBuffer(&plist->imask);
366 freeDeviceBuffer(&plist->excl);
368 if (nb->bUseTwoStreams)
370 auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
371 freeDeviceBuffer(&plist_nl->sci);
372 freeDeviceBuffer(&plist_nl->cj4);
373 freeDeviceBuffer(&plist_nl->imask);
374 freeDeviceBuffer(&plist_nl->excl);
380 nb->nbst.eLJ = nullptr;
382 pfree(nb->nbst.eElec);
383 nb->nbst.eElec = nullptr;
385 pfree(nb->nbst.fShift);
386 nb->nbst.fShift = nullptr;
395 fprintf(debug, "Cleaned up CUDA data structures.\n");
399 int gpu_min_ci_balanced(NbnxmGpu* nb)
401 return nb != nullptr ? gpu_min_ci_balanced_factor * nb->deviceContext_->deviceInfo().prop.multiProcessorCount
405 void* gpu_get_xq(NbnxmGpu* nb)
409 return static_cast<void*>(nb->atdat->xq);
412 DeviceBuffer<gmx::RVec> gpu_get_f(NbnxmGpu* nb)
416 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->f);
419 DeviceBuffer<gmx::RVec> gpu_get_fshift(NbnxmGpu* nb)
423 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->fShift);
426 /* Initialization for X buffer operations on GPU. */
427 /* TODO Remove explicit pinning from host arrays from here and manage in a more natural way*/
428 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet& gridSet, NbnxmGpu* gpu_nbv)
430 const DeviceStream& localStream = *gpu_nbv->deviceStreams[InteractionLocality::Local];
431 bool bDoTime = gpu_nbv->bDoTime;
432 const int maxNumColumns = gridSet.numColumnsMax();
434 reallocateDeviceBuffer(&gpu_nbv->cxy_na,
435 maxNumColumns * gridSet.grids().size(),
437 &gpu_nbv->ncxy_na_alloc,
438 *gpu_nbv->deviceContext_);
439 reallocateDeviceBuffer(&gpu_nbv->cxy_ind,
440 maxNumColumns * gridSet.grids().size(),
442 &gpu_nbv->ncxy_ind_alloc,
443 *gpu_nbv->deviceContext_);
445 for (unsigned int g = 0; g < gridSet.grids().size(); g++)
448 const Nbnxm::Grid& grid = gridSet.grids()[g];
450 const int numColumns = grid.numColumns();
451 const int* atomIndices = gridSet.atomIndices().data();
452 const int atomIndicesSize = gridSet.atomIndices().size();
453 const int* cxy_na = grid.cxy_na().data();
454 const int* cxy_ind = grid.cxy_ind().data();
456 reallocateDeviceBuffer(&gpu_nbv->atomIndices,
458 &gpu_nbv->atomIndicesSize,
459 &gpu_nbv->atomIndicesSize_alloc,
460 *gpu_nbv->deviceContext_);
462 if (atomIndicesSize > 0)
467 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(localStream);
470 copyToDeviceBuffer(&gpu_nbv->atomIndices,
475 GpuApiCallBehavior::Async,
480 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(localStream);
488 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(localStream);
491 int* destPtr = &gpu_nbv->cxy_na[maxNumColumns * g];
493 &destPtr, cxy_na, 0, numColumns, localStream, GpuApiCallBehavior::Async, nullptr);
497 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(localStream);
502 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(localStream);
505 destPtr = &gpu_nbv->cxy_ind[maxNumColumns * g];
507 &destPtr, cxy_ind, 0, numColumns, localStream, GpuApiCallBehavior::Async, nullptr);
511 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(localStream);
516 if (gpu_nbv->bUseTwoStreams)
518 // The above data is transferred on the local stream but is a
519 // dependency of the nonlocal stream (specifically the nonlocal X
520 // buf ops kernel). We therefore set a dependency to ensure
521 // that the nonlocal stream waits on the local stream here.
522 // This call records an event in the local stream:
523 gpu_nbv->misc_ops_and_local_H2D_done.markEvent(
524 *gpu_nbv->deviceStreams[Nbnxm::InteractionLocality::Local]);
525 // ...and this call instructs the nonlocal stream to wait on that event:
526 gpu_nbv->misc_ops_and_local_H2D_done.enqueueWaitEvent(
527 *gpu_nbv->deviceStreams[Nbnxm::InteractionLocality::NonLocal]);