Improve workload data structures' docs
[gromacs.git] / src / gromacs / listed_forces / gpubonded_impl.cu
blobe71738f586c160bbb41ef7896b9eadff3812049c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief Implements GPU bonded lists for CUDA
38  *
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>
43  *
44  * \ingroup module_listed_forces
45  */
47 #include "gmxpre.h"
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/devicebuffer.h"
54 #include "gromacs/mdtypes/enerdata.h"
55 #include "gromacs/timing/wallcycle.h"
56 #include "gromacs/topology/forcefieldparameters.h"
58 struct t_forcerec;
60 namespace gmx
63 // ---- GpuBonded::Impl
65 GpuBonded::Impl::Impl(const gmx_ffparams_t &ffparams,
66                       void                 *streamPtr,
67                       gmx_wallcycle        *wcycle)
69     stream_ = *static_cast<CommandStream*>(streamPtr);
70     wcycle_ = wcycle;
72     allocateDeviceBuffer(&d_forceParams_, ffparams.numTypes(), nullptr);
73     // This could be an async transfer (if the source is pinned), so
74     // long as it uses the same stream as the kernels and we are happy
75     // to consume additional pinned pages.
76     copyToDeviceBuffer(&d_forceParams_, ffparams.iparams.data(),
77                        0, ffparams.numTypes(),
78                        stream_, GpuApiCallBehavior::Sync, nullptr);
79     vTot_.resize(F_NRE);
80     allocateDeviceBuffer(&d_vTot_, F_NRE, nullptr);
81     clearDeviceBufferAsync(&d_vTot_, 0, F_NRE, stream_);
83     for (int fType = 0; fType < F_NRE; fType++)
84     {
85         d_iLists_[fType].nr     = 0;
86         d_iLists_[fType].iatoms = nullptr;
87         d_iLists_[fType].nalloc = 0;
88     }
90     kernelParams_.d_forceParams          = d_forceParams_;
91     kernelParams_.d_xq                   = d_xq_;
92     kernelParams_.d_f                    = d_f_;
93     kernelParams_.d_fShift               = d_fShift_;
94     kernelParams_.d_vTot                 = d_vTot_;
95     for (int i = 0; i < numFTypesOnGpu; i++)
96     {
97         kernelParams_.d_iatoms[i]        = nullptr;
98         kernelParams_.fTypeRangeStart[i] = 0;
99         kernelParams_.fTypeRangeEnd[i]   = -1;
100     }
103 GpuBonded::Impl::~Impl()
105     for (int fType : fTypesOnGpu)
106     {
107         if (d_iLists_[fType].iatoms)
108         {
109             freeDeviceBuffer(&d_iLists_[fType].iatoms);
110             d_iLists_[fType].iatoms = nullptr;
111         }
112     }
114     freeDeviceBuffer(&d_forceParams_);
115     freeDeviceBuffer(&d_vTot_);
118 //! Return whether function type \p fType in \p idef has perturbed interactions
119 static bool fTypeHasPerturbedEntries(const t_idef  &idef,
120                                      int            fType)
122     GMX_ASSERT(idef.ilsort == ilsortNO_FE || idef.ilsort == ilsortFE_SORTED,
123                "Perturbed interations should be sorted here");
125     const t_ilist &ilist = idef.il[fType];
127     return (idef.ilsort != ilsortNO_FE && ilist.nr_nonperturbed != ilist.nr);
130 //! Converts \p src with atom indices in state order to \p dest in nbnxn order
131 static void convertIlistToNbnxnOrder(const t_ilist       &src,
132                                      HostInteractionList *dest,
133                                      int                  numAtomsPerInteraction,
134                                      ArrayRef<const int>  nbnxnAtomOrder)
136     GMX_ASSERT(src.size() == 0 || !nbnxnAtomOrder.empty(), "We need the nbnxn atom order");
138     dest->iatoms.resize(src.size());
140     // TODO use OpenMP to parallelise this loop
141     for (int i = 0; i < src.size(); i += 1 + numAtomsPerInteraction)
142     {
143         dest->iatoms[i] = src.iatoms[i];
144         for (int a = 0; a < numAtomsPerInteraction; a++)
145         {
146             dest->iatoms[i + 1 + a] = nbnxnAtomOrder[src.iatoms[i + 1 + a]];
147         }
148     }
151 //! Returns \p input rounded up to the closest multiple of \p factor.
152 static inline int roundUpToFactor(const int input, const int factor)
154     GMX_ASSERT(factor > 0, "The factor to round up to must be > 0.");
156     int remainder = input % factor;
158     if (remainder == 0)
159     {
160         return (input);
161     }
162     return (input + (factor - remainder));
165 // TODO Consider whether this function should be a factory method that
166 // makes an object that is the only one capable of the device
167 // operations needed for the lifetime of an interaction list. This
168 // would be harder to misuse than GpuBonded, and exchange the problem
169 // of naming this method for the problem of what to name the
170 // BondedDeviceInteractionListHandler type.
172 /*! Divides bonded interactions over threads and GPU.
173  *  The bonded interactions are assigned by interaction type to GPU threads. The intereaction
174  *  types are assigned in blocks sized as <warp_size>. The beginning and end (thread index) of each
175  *  interaction type are stored in kernelParams_. Pointers to the relevant data structures on the
176  *  GPU are also stored in kernelParams_.
177  */
178 void
179 GpuBonded::Impl::updateInteractionListsAndDeviceBuffers(ArrayRef<const int>  nbnxnAtomOrder,
180                                                         const t_idef        &idef,
181                                                         void                *d_xqPtr,
182                                                         void                *d_fPtr,
183                                                         void                *d_fShiftPtr)
185     // TODO wallcycle sub start
186     haveInteractions_ = false;
187     int fTypesCounter = 0;
189     for (int fType : fTypesOnGpu)
190     {
191         auto &iList = iLists_[fType];
193         /* Perturbation is not implemented in the GPU bonded kernels.
194          * But instead of doing all interactions on the CPU, we can
195          * still easily handle the types that have no perturbed
196          * interactions on the GPU. */
197         if (idef.il[fType].nr > 0 && !fTypeHasPerturbedEntries(idef, fType))
198         {
199             haveInteractions_ = true;
201             convertIlistToNbnxnOrder(idef.il[fType],
202                                      &iList,
203                                      NRAL(fType), nbnxnAtomOrder);
204         }
205         else
206         {
207             iList.iatoms.clear();
208         }
210         // Update the device if necessary. This can leave some
211         // allocations on the device when the host size decreases to
212         // zero, which is OK, since we deallocate everything at the
213         // end.
214         if (iList.size() > 0)
215         {
216             t_ilist &d_iList = d_iLists_[fType];
218             reallocateDeviceBuffer(&d_iList.iatoms, iList.size(), &d_iList.nr, &d_iList.nalloc, nullptr);
220             copyToDeviceBuffer(&d_iList.iatoms, iList.iatoms.data(),
221                                0, iList.size(),
222                                stream_, GpuApiCallBehavior::Async, nullptr);
223         }
224         kernelParams_.fTypesOnGpu[fTypesCounter]    = fType;
225         kernelParams_.numFTypeIAtoms[fTypesCounter] = iList.size();
226         int numBonds = iList.size() / (interaction_function[fType].nratoms + 1);
227         kernelParams_.numFTypeBonds[fTypesCounter]  = numBonds;
228         kernelParams_.d_iatoms[fTypesCounter]       = d_iLists_[fType].iatoms;
229         if (fTypesCounter == 0)
230         {
231             kernelParams_.fTypeRangeStart[fTypesCounter] = 0;
232         }
233         else
234         {
235             kernelParams_.fTypeRangeStart[fTypesCounter] = kernelParams_.fTypeRangeEnd[fTypesCounter - 1] + 1;
236         }
237         kernelParams_.fTypeRangeEnd[fTypesCounter] = kernelParams_.fTypeRangeStart[fTypesCounter] + roundUpToFactor(numBonds, warp_size) - 1;
239         GMX_ASSERT(numBonds > 0 || kernelParams_.fTypeRangeEnd[fTypesCounter] <= kernelParams_.fTypeRangeStart[fTypesCounter],
240                    "Invalid GPU listed forces setup. numBonds must be > 0 if there are threads allocated to do work on that interaction function type.");
241         GMX_ASSERT(kernelParams_.fTypeRangeStart[fTypesCounter] % warp_size == 0 && (kernelParams_.fTypeRangeEnd[fTypesCounter] + 1) % warp_size == 0,
242                    "The bonded interactions must be assigned to the GPU in blocks of warp size.");
244         fTypesCounter++;
245     }
247     d_xq_                           = static_cast<float4 *>(d_xqPtr);
248     d_f_                            = static_cast<fvec *>(d_fPtr);
249     d_fShift_                       = static_cast<fvec *>(d_fShiftPtr);
251     kernelParams_.d_xq              = d_xq_;
252     kernelParams_.d_f               = d_f_;
253     kernelParams_.d_fShift          = d_fShift_;
254     kernelParams_.d_forceParams     = d_forceParams_;
255     kernelParams_.d_vTot            = d_vTot_;
257     // TODO wallcycle sub stop
260 bool
261 GpuBonded::Impl::haveInteractions() const
263     return haveInteractions_;
266 void
267 GpuBonded::Impl::launchEnergyTransfer()
269     GMX_ASSERT(haveInteractions_, "No GPU bonded interactions, so no energies will be computed, so transfer should not be called");
271     wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_GPU_BONDED);
272     // TODO add conditional on whether there has been any compute (and make sure host buffer doesn't contain garbage)
273     float *h_vTot   = vTot_.data();
274     copyFromDeviceBuffer(h_vTot, &d_vTot_,
275                          0, F_NRE,
276                          stream_, GpuApiCallBehavior::Async, nullptr);
277     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_BONDED);
280 void
281 GpuBonded::Impl::waitAccumulateEnergyTerms(gmx_enerdata_t *enerd)
283     GMX_ASSERT(haveInteractions_, "No GPU bonded interactions, so no energies will be computed or transferred, so accumulation should not occur");
285     wallcycle_start(wcycle_, ewcWAIT_GPU_BONDED);
286     cudaError_t stat = cudaStreamSynchronize(stream_);
287     CU_RET_ERR(stat, "D2H transfer of bonded energies failed");
288     wallcycle_stop(wcycle_, ewcWAIT_GPU_BONDED);
290     for (int fType : fTypesOnGpu)
291     {
292         if (fType != F_LJ14 && fType != F_COUL14)
293         {
294             enerd->term[fType] += vTot_[fType];
295         }
296     }
298     // Note: We do not support energy groups here
299     gmx_grppairener_t *grppener = &enerd->grpp;
300     GMX_RELEASE_ASSERT(grppener->nener == 1, "No energy group support for bondeds on the GPU");
301     grppener->ener[egLJ14][0]   += vTot_[F_LJ14];
302     grppener->ener[egCOUL14][0] += vTot_[F_COUL14];
305 void
306 GpuBonded::Impl::clearEnergies()
308     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
309     wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_GPU_BONDED);
310     clearDeviceBufferAsync(&d_vTot_, 0, F_NRE, stream_);
311     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_BONDED);
312     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
315 // ---- GpuBonded
317 GpuBonded::GpuBonded(const gmx_ffparams_t &ffparams,
318                      void                 *streamPtr,
319                      gmx_wallcycle        *wcycle)
320     : impl_(new Impl(ffparams, streamPtr, wcycle))
324 GpuBonded::~GpuBonded() = default;
326 void
327 GpuBonded::updateInteractionListsAndDeviceBuffers(ArrayRef<const int>  nbnxnAtomOrder,
328                                                   const t_idef        &idef,
329                                                   void                *d_xq,
330                                                   void                *d_f,
331                                                   void                *d_fShift)
333     impl_->updateInteractionListsAndDeviceBuffers
334         (nbnxnAtomOrder, idef, d_xq, d_f, d_fShift);
337 bool
338 GpuBonded::haveInteractions() const
340     return impl_->haveInteractions();
343 void
344 GpuBonded::launchEnergyTransfer()
346     impl_->launchEnergyTransfer();
349 void
350 GpuBonded::waitAccumulateEnergyTerms(gmx_enerdata_t *enerd)
352     impl_->waitAccumulateEnergyTerms(enerd);
355 void
356 GpuBonded::clearEnergies()
358     impl_->clearEnergies();
361 }   // namespace gmx