2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * \brief Implements GPU bonded lists for CUDA
39 * \author Berk Hess <hess@kth.se>
40 * \author Szilárd Páll <pall.szilard@gmail.com>
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \author Magnus Lundborg <lundborg.magnus@gmail.com>
44 * \ingroup module_listed_forces
49 #include "gpubonded_impl.h"
51 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
52 #include "gromacs/gpu_utils/cudautils.cuh"
53 #include "gromacs/gpu_utils/device_context.h"
54 #include "gromacs/gpu_utils/device_stream.h"
55 #include "gromacs/gpu_utils/devicebuffer.h"
56 #include "gromacs/gpu_utils/typecasts.cuh"
57 #include "gromacs/mdtypes/enerdata.h"
58 #include "gromacs/timing/wallcycle.h"
59 #include "gromacs/topology/forcefieldparameters.h"
65 // Number of CUDA threads in a block
66 constexpr static int c_threadsPerBlock = 256;
68 // ---- GpuBonded::Impl
70 GpuBonded::Impl::Impl(const gmx_ffparams_t& ffparams,
71 const float electrostaticsScaleFactor,
72 const DeviceContext& deviceContext,
73 const DeviceStream& deviceStream,
74 gmx_wallcycle* wcycle) :
75 deviceContext_(deviceContext),
76 deviceStream_(deviceStream)
78 GMX_RELEASE_ASSERT(deviceStream.isValid(),
79 "Can't run GPU version of bonded forces in stream that is not valid.");
82 c_threadsPerBlock >= c_numShiftVectors,
83 "Threads per block in GPU bonded must be >= c_numShiftVectors for the virial kernel "
88 allocateDeviceBuffer(&d_forceParams_, ffparams.numTypes(), deviceContext_);
89 // This could be an async transfer (if the source is pinned), so
90 // long as it uses the same stream as the kernels and we are happy
91 // to consume additional pinned pages.
92 copyToDeviceBuffer(&d_forceParams_,
93 ffparams.iparams.data(),
97 GpuApiCallBehavior::Sync,
100 allocateDeviceBuffer(&d_vTot_, F_NRE, deviceContext_);
101 clearDeviceBufferAsync(&d_vTot_, 0, F_NRE, deviceStream_);
103 kernelParams_.electrostaticsScaleFactor = electrostaticsScaleFactor;
104 kernelParams_.d_forceParams = d_forceParams_;
105 kernelParams_.d_xq = d_xq_;
106 kernelParams_.d_f = d_f_;
107 kernelParams_.d_fShift = d_fShift_;
108 kernelParams_.d_vTot = d_vTot_;
109 for (int i = 0; i < numFTypesOnGpu; i++)
111 kernelParams_.d_iatoms[i] = nullptr;
112 kernelParams_.fTypeRangeStart[i] = 0;
113 kernelParams_.fTypeRangeEnd[i] = -1;
116 int fTypeRangeEnd = kernelParams_.fTypeRangeEnd[numFTypesOnGpu - 1];
118 kernelLaunchConfig_.blockSize[0] = c_threadsPerBlock;
119 kernelLaunchConfig_.blockSize[1] = 1;
120 kernelLaunchConfig_.blockSize[2] = 1;
121 kernelLaunchConfig_.gridSize[0] = (fTypeRangeEnd + c_threadsPerBlock) / c_threadsPerBlock;
122 kernelLaunchConfig_.gridSize[1] = 1;
123 kernelLaunchConfig_.gridSize[2] = 1;
124 kernelLaunchConfig_.sharedMemorySize =
125 c_numShiftVectors * sizeof(float3) + (c_threadsPerBlock / warp_size) * 3 * sizeof(float);
128 GpuBonded::Impl::~Impl()
130 for (int fType : fTypesOnGpu)
132 if (d_iLists_[fType].iatoms)
134 freeDeviceBuffer(&d_iLists_[fType].iatoms);
135 d_iLists_[fType].iatoms = nullptr;
139 freeDeviceBuffer(&d_forceParams_);
140 freeDeviceBuffer(&d_vTot_);
143 //! Return whether function type \p fType in \p idef has perturbed interactions
144 static bool fTypeHasPerturbedEntries(const InteractionDefinitions& idef, int fType)
146 GMX_ASSERT(idef.ilsort == ilsortNO_FE || idef.ilsort == ilsortFE_SORTED,
147 "Perturbed interations should be sorted here");
149 const InteractionList& ilist = idef.il[fType];
151 return (idef.ilsort != ilsortNO_FE && idef.numNonperturbedInteractions[fType] != ilist.size());
154 //! Converts \p src with atom indices in state order to \p dest in nbnxn order
155 static void convertIlistToNbnxnOrder(const InteractionList& src,
156 HostInteractionList* dest,
157 int numAtomsPerInteraction,
158 ArrayRef<const int> nbnxnAtomOrder)
160 GMX_ASSERT(src.size() == 0 || !nbnxnAtomOrder.empty(), "We need the nbnxn atom order");
162 dest->iatoms.resize(src.size());
164 // TODO use OpenMP to parallelise this loop
165 for (int i = 0; i < src.size(); i += 1 + numAtomsPerInteraction)
167 dest->iatoms[i] = src.iatoms[i];
168 for (int a = 0; a < numAtomsPerInteraction; a++)
170 dest->iatoms[i + 1 + a] = nbnxnAtomOrder[src.iatoms[i + 1 + a]];
175 //! Returns \p input rounded up to the closest multiple of \p factor.
176 static inline int roundUpToFactor(const int input, const int factor)
178 GMX_ASSERT(factor > 0, "The factor to round up to must be > 0.");
180 int remainder = input % factor;
186 return (input + (factor - remainder));
189 // TODO Consider whether this function should be a factory method that
190 // makes an object that is the only one capable of the device
191 // operations needed for the lifetime of an interaction list. This
192 // would be harder to misuse than GpuBonded, and exchange the problem
193 // of naming this method for the problem of what to name the
194 // BondedDeviceInteractionListHandler type.
196 /*! Divides bonded interactions over threads and GPU.
197 * The bonded interactions are assigned by interaction type to GPU threads. The intereaction
198 * types are assigned in blocks sized as <warp_size>. The beginning and end (thread index) of each
199 * interaction type are stored in kernelParams_. Pointers to the relevant data structures on the
200 * GPU are also stored in kernelParams_.
202 * \todo Use DeviceBuffer for the d_xqPtr.
204 void GpuBonded::Impl::updateInteractionListsAndDeviceBuffers(ArrayRef<const int> nbnxnAtomOrder,
205 const InteractionDefinitions& idef,
207 DeviceBuffer<RVec> d_fPtr,
208 DeviceBuffer<RVec> d_fShiftPtr)
210 // TODO wallcycle sub start
211 haveInteractions_ = false;
212 int fTypesCounter = 0;
214 for (int fType : fTypesOnGpu)
216 auto& iList = iLists_[fType];
218 /* Perturbation is not implemented in the GPU bonded kernels.
219 * But instead of doing all interactions on the CPU, we can
220 * still easily handle the types that have no perturbed
221 * interactions on the GPU. */
222 if (!idef.il[fType].empty() && !fTypeHasPerturbedEntries(idef, fType))
224 haveInteractions_ = true;
226 convertIlistToNbnxnOrder(idef.il[fType], &iList, NRAL(fType), nbnxnAtomOrder);
230 iList.iatoms.clear();
233 // Update the device if necessary. This can leave some
234 // allocations on the device when the host size decreases to
235 // zero, which is OK, since we deallocate everything at the
237 if (iList.size() > 0)
239 t_ilist& d_iList = d_iLists_[fType];
241 reallocateDeviceBuffer(
242 &d_iList.iatoms, iList.size(), &d_iList.nr, &d_iList.nalloc, deviceContext_);
244 copyToDeviceBuffer(&d_iList.iatoms,
249 GpuApiCallBehavior::Async,
252 kernelParams_.fTypesOnGpu[fTypesCounter] = fType;
253 kernelParams_.numFTypeIAtoms[fTypesCounter] = iList.size();
254 int numBonds = iList.size() / (interaction_function[fType].nratoms + 1);
255 kernelParams_.numFTypeBonds[fTypesCounter] = numBonds;
256 kernelParams_.d_iatoms[fTypesCounter] = d_iLists_[fType].iatoms;
257 if (fTypesCounter == 0)
259 kernelParams_.fTypeRangeStart[fTypesCounter] = 0;
263 kernelParams_.fTypeRangeStart[fTypesCounter] =
264 kernelParams_.fTypeRangeEnd[fTypesCounter - 1] + 1;
266 kernelParams_.fTypeRangeEnd[fTypesCounter] = kernelParams_.fTypeRangeStart[fTypesCounter]
267 + roundUpToFactor(numBonds, warp_size) - 1;
269 GMX_ASSERT(numBonds > 0
270 || kernelParams_.fTypeRangeEnd[fTypesCounter]
271 <= kernelParams_.fTypeRangeStart[fTypesCounter],
272 "Invalid GPU listed forces setup. numBonds must be > 0 if there are threads "
273 "allocated to do work on that interaction function type.");
274 GMX_ASSERT(kernelParams_.fTypeRangeStart[fTypesCounter] % warp_size == 0
275 && (kernelParams_.fTypeRangeEnd[fTypesCounter] + 1) % warp_size == 0,
276 "The bonded interactions must be assigned to the GPU in blocks of warp size.");
281 int fTypeRangeEnd = kernelParams_.fTypeRangeEnd[numFTypesOnGpu - 1];
282 kernelLaunchConfig_.gridSize[0] = (fTypeRangeEnd + c_threadsPerBlock) / c_threadsPerBlock;
284 d_xq_ = static_cast<float4*>(d_xqPtr);
285 d_f_ = asFloat3(d_fPtr);
286 d_fShift_ = asFloat3(d_fShiftPtr);
288 kernelParams_.d_xq = d_xq_;
289 kernelParams_.d_f = d_f_;
290 kernelParams_.d_fShift = d_fShift_;
291 kernelParams_.d_forceParams = d_forceParams_;
292 kernelParams_.d_vTot = d_vTot_;
294 // TODO wallcycle sub stop
297 void GpuBonded::Impl::setPbc(PbcType pbcType, const matrix box, bool canMoleculeSpanPbc)
300 setPbcAiuc(canMoleculeSpanPbc ? numPbcDimensions(pbcType) : 0, box, &pbcAiuc);
301 kernelParams_.pbcAiuc = pbcAiuc;
304 bool GpuBonded::Impl::haveInteractions() const
306 return haveInteractions_;
309 void GpuBonded::Impl::launchEnergyTransfer()
311 GMX_ASSERT(haveInteractions_,
312 "No GPU bonded interactions, so no energies will be computed, so transfer should "
314 wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_GPU_BONDED);
315 // TODO add conditional on whether there has been any compute (and make sure host buffer doesn't contain garbage)
316 float* h_vTot = vTot_.data();
317 copyFromDeviceBuffer(h_vTot, &d_vTot_, 0, F_NRE, deviceStream_, GpuApiCallBehavior::Async, nullptr);
318 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_BONDED);
321 void GpuBonded::Impl::waitAccumulateEnergyTerms(gmx_enerdata_t* enerd)
323 GMX_ASSERT(haveInteractions_,
324 "No GPU bonded interactions, so no energies will be computed or transferred, so "
325 "accumulation should not occur");
327 wallcycle_start(wcycle_, ewcWAIT_GPU_BONDED);
328 cudaError_t stat = cudaStreamSynchronize(deviceStream_.stream());
329 CU_RET_ERR(stat, "D2H transfer of bonded energies failed");
330 wallcycle_stop(wcycle_, ewcWAIT_GPU_BONDED);
332 for (int fType : fTypesOnGpu)
334 if (fType != F_LJ14 && fType != F_COUL14)
336 enerd->term[fType] += vTot_[fType];
340 // Note: We do not support energy groups here
341 gmx_grppairener_t* grppener = &enerd->grpp;
342 GMX_RELEASE_ASSERT(grppener->nener == 1, "No energy group support for bondeds on the GPU");
343 grppener->energyGroupPairTerms[NonBondedEnergyTerms::LJ14][0] += vTot_[F_LJ14];
344 grppener->energyGroupPairTerms[NonBondedEnergyTerms::Coulomb14][0] += vTot_[F_COUL14];
347 void GpuBonded::Impl::clearEnergies()
349 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
350 wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_GPU_BONDED);
351 clearDeviceBufferAsync(&d_vTot_, 0, F_NRE, deviceStream_);
352 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_BONDED);
353 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
358 GpuBonded::GpuBonded(const gmx_ffparams_t& ffparams,
359 const float electrostaticsScaleFactor,
360 const DeviceContext& deviceContext,
361 const DeviceStream& deviceStream,
362 gmx_wallcycle* wcycle) :
363 impl_(new Impl(ffparams, electrostaticsScaleFactor, deviceContext, deviceStream, wcycle))
367 GpuBonded::~GpuBonded() = default;
369 void GpuBonded::updateInteractionListsAndDeviceBuffers(ArrayRef<const int> nbnxnAtomOrder,
370 const InteractionDefinitions& idef,
372 DeviceBuffer<RVec> d_f,
373 DeviceBuffer<RVec> d_fShift)
375 impl_->updateInteractionListsAndDeviceBuffers(nbnxnAtomOrder, idef, d_xq, d_f, d_fShift);
378 void GpuBonded::setPbc(PbcType pbcType, const matrix box, bool canMoleculeSpanPbc)
380 impl_->setPbc(pbcType, box, canMoleculeSpanPbc);
383 bool GpuBonded::haveInteractions() const
385 return impl_->haveInteractions();
388 void GpuBonded::setPbcAndlaunchKernel(PbcType pbcType,
390 bool canMoleculeSpanPbc,
391 const gmx::StepWorkload& stepWork)
393 setPbc(pbcType, box, canMoleculeSpanPbc);
394 launchKernel(stepWork);
397 void GpuBonded::launchEnergyTransfer()
399 impl_->launchEnergyTransfer();
402 void GpuBonded::waitAccumulateEnergyTerms(gmx_enerdata_t* enerd)
404 impl_->waitAccumulateEnergyTerms(enerd);
407 void GpuBonded::clearEnergies()
409 impl_->clearEnergies();