Make use of the DeviceStreamManager
[alexxy/gromacs.git] / src / gromacs / listed_forces / gpubonded_impl.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020, 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  */
46
47 #include "gmxpre.h"
48
49 #include "gpubonded_impl.h"
50
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"
60
61 struct t_forcerec;
62
63 namespace gmx
64 {
65
66 // ---- GpuBonded::Impl
67
68 GpuBonded::Impl::Impl(const gmx_ffparams_t& ffparams,
69                       const DeviceContext&  deviceContext,
70                       const DeviceStream&   deviceStream,
71                       gmx_wallcycle*        wcycle) :
72     deviceContext_(deviceContext),
73     deviceStream_(deviceStream)
74 {
75     GMX_RELEASE_ASSERT(deviceStream.isValid(),
76                        "Can't run GPU version of bonded forces in stream that is not valid.");
77
78     wcycle_ = wcycle;
79
80     allocateDeviceBuffer(&d_forceParams_, ffparams.numTypes(), deviceContext_);
81     // This could be an async transfer (if the source is pinned), so
82     // long as it uses the same stream as the kernels and we are happy
83     // to consume additional pinned pages.
84     copyToDeviceBuffer(&d_forceParams_, ffparams.iparams.data(), 0, ffparams.numTypes(),
85                        deviceStream_, GpuApiCallBehavior::Sync, nullptr);
86     vTot_.resize(F_NRE);
87     allocateDeviceBuffer(&d_vTot_, F_NRE, deviceContext_);
88     clearDeviceBufferAsync(&d_vTot_, 0, F_NRE, deviceStream_);
89
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     }
101 }
102
103 GpuBonded::Impl::~Impl()
104 {
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     }
113
114     freeDeviceBuffer(&d_forceParams_);
115     freeDeviceBuffer(&d_vTot_);
116 }
117
118 //! Return whether function type \p fType in \p idef has perturbed interactions
119 static bool fTypeHasPerturbedEntries(const InteractionDefinitions& idef, int fType)
120 {
121     GMX_ASSERT(idef.ilsort == ilsortNO_FE || idef.ilsort == ilsortFE_SORTED,
122                "Perturbed interations should be sorted here");
123
124     const InteractionList& ilist = idef.il[fType];
125
126     return (idef.ilsort != ilsortNO_FE && idef.numNonperturbedInteractions[fType] != ilist.size());
127 }
128
129 //! Converts \p src with atom indices in state order to \p dest in nbnxn order
130 static void convertIlistToNbnxnOrder(const InteractionList& src,
131                                      HostInteractionList*   dest,
132                                      int                    numAtomsPerInteraction,
133                                      ArrayRef<const int>    nbnxnAtomOrder)
134 {
135     GMX_ASSERT(src.size() == 0 || !nbnxnAtomOrder.empty(), "We need the nbnxn atom order");
136
137     dest->iatoms.resize(src.size());
138
139     // TODO use OpenMP to parallelise this loop
140     for (int i = 0; i < src.size(); i += 1 + numAtomsPerInteraction)
141     {
142         dest->iatoms[i] = src.iatoms[i];
143         for (int a = 0; a < numAtomsPerInteraction; a++)
144         {
145             dest->iatoms[i + 1 + a] = nbnxnAtomOrder[src.iatoms[i + 1 + a]];
146         }
147     }
148 }
149
150 //! Returns \p input rounded up to the closest multiple of \p factor.
151 static inline int roundUpToFactor(const int input, const int factor)
152 {
153     GMX_ASSERT(factor > 0, "The factor to round up to must be > 0.");
154
155     int remainder = input % factor;
156
157     if (remainder == 0)
158     {
159         return (input);
160     }
161     return (input + (factor - remainder));
162 }
163
164 // TODO Consider whether this function should be a factory method that
165 // makes an object that is the only one capable of the device
166 // operations needed for the lifetime of an interaction list. This
167 // would be harder to misuse than GpuBonded, and exchange the problem
168 // of naming this method for the problem of what to name the
169 // BondedDeviceInteractionListHandler type.
170
171 /*! Divides bonded interactions over threads and GPU.
172  *  The bonded interactions are assigned by interaction type to GPU threads. The intereaction
173  *  types are assigned in blocks sized as <warp_size>. The beginning and end (thread index) of each
174  *  interaction type are stored in kernelParams_. Pointers to the relevant data structures on the
175  *  GPU are also stored in kernelParams_.
176  *
177  * \todo Use DeviceBuffer for the d_xqPtr.
178  */
179 void GpuBonded::Impl::updateInteractionListsAndDeviceBuffers(ArrayRef<const int> nbnxnAtomOrder,
180                                                              const InteractionDefinitions& idef,
181                                                              void*                         d_xqPtr,
182                                                              DeviceBuffer<RVec>            d_fPtr,
183                                                              DeviceBuffer<RVec> d_fShiftPtr)
184 {
185     // TODO wallcycle sub start
186     haveInteractions_ = false;
187     int fTypesCounter = 0;
188
189     for (int fType : fTypesOnGpu)
190     {
191         auto& iList = iLists_[fType];
192
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].empty() && !fTypeHasPerturbedEntries(idef, fType))
198         {
199             haveInteractions_ = true;
200
201             convertIlistToNbnxnOrder(idef.il[fType], &iList, NRAL(fType), nbnxnAtomOrder);
202         }
203         else
204         {
205             iList.iatoms.clear();
206         }
207
208         // Update the device if necessary. This can leave some
209         // allocations on the device when the host size decreases to
210         // zero, which is OK, since we deallocate everything at the
211         // end.
212         if (iList.size() > 0)
213         {
214             t_ilist& d_iList = d_iLists_[fType];
215
216             reallocateDeviceBuffer(&d_iList.iatoms, iList.size(), &d_iList.nr, &d_iList.nalloc,
217                                    deviceContext_);
218
219             copyToDeviceBuffer(&d_iList.iatoms, iList.iatoms.data(), 0, iList.size(), deviceStream_,
220                                GpuApiCallBehavior::Async, nullptr);
221         }
222         kernelParams_.fTypesOnGpu[fTypesCounter]    = fType;
223         kernelParams_.numFTypeIAtoms[fTypesCounter] = iList.size();
224         int numBonds = iList.size() / (interaction_function[fType].nratoms + 1);
225         kernelParams_.numFTypeBonds[fTypesCounter] = numBonds;
226         kernelParams_.d_iatoms[fTypesCounter]      = d_iLists_[fType].iatoms;
227         if (fTypesCounter == 0)
228         {
229             kernelParams_.fTypeRangeStart[fTypesCounter] = 0;
230         }
231         else
232         {
233             kernelParams_.fTypeRangeStart[fTypesCounter] =
234                     kernelParams_.fTypeRangeEnd[fTypesCounter - 1] + 1;
235         }
236         kernelParams_.fTypeRangeEnd[fTypesCounter] = kernelParams_.fTypeRangeStart[fTypesCounter]
237                                                      + roundUpToFactor(numBonds, warp_size) - 1;
238
239         GMX_ASSERT(numBonds > 0
240                            || kernelParams_.fTypeRangeEnd[fTypesCounter]
241                                       <= kernelParams_.fTypeRangeStart[fTypesCounter],
242                    "Invalid GPU listed forces setup. numBonds must be > 0 if there are threads "
243                    "allocated to do work on that interaction function type.");
244         GMX_ASSERT(kernelParams_.fTypeRangeStart[fTypesCounter] % warp_size == 0
245                            && (kernelParams_.fTypeRangeEnd[fTypesCounter] + 1) % warp_size == 0,
246                    "The bonded interactions must be assigned to the GPU in blocks of warp size.");
247
248         fTypesCounter++;
249     }
250
251     d_xq_     = static_cast<float4*>(d_xqPtr);
252     d_f_      = asFloat3(d_fPtr);
253     d_fShift_ = asFloat3(d_fShiftPtr);
254
255     kernelParams_.d_xq          = d_xq_;
256     kernelParams_.d_f           = d_f_;
257     kernelParams_.d_fShift      = d_fShift_;
258     kernelParams_.d_forceParams = d_forceParams_;
259     kernelParams_.d_vTot        = d_vTot_;
260
261     // TODO wallcycle sub stop
262 }
263
264 bool GpuBonded::Impl::haveInteractions() const
265 {
266     return haveInteractions_;
267 }
268
269 void GpuBonded::Impl::launchEnergyTransfer()
270 {
271     GMX_ASSERT(haveInteractions_,
272                "No GPU bonded interactions, so no energies will be computed, so transfer should "
273                "not be called");
274     wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_GPU_BONDED);
275     // TODO add conditional on whether there has been any compute (and make sure host buffer doesn't contain garbage)
276     float* h_vTot = vTot_.data();
277     copyFromDeviceBuffer(h_vTot, &d_vTot_, 0, F_NRE, deviceStream_, GpuApiCallBehavior::Async, nullptr);
278     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_BONDED);
279 }
280
281 void GpuBonded::Impl::waitAccumulateEnergyTerms(gmx_enerdata_t* enerd)
282 {
283     GMX_ASSERT(haveInteractions_,
284                "No GPU bonded interactions, so no energies will be computed or transferred, so "
285                "accumulation should not occur");
286
287     wallcycle_start(wcycle_, ewcWAIT_GPU_BONDED);
288     cudaError_t stat = cudaStreamSynchronize(deviceStream_.stream());
289     CU_RET_ERR(stat, "D2H transfer of bonded energies failed");
290     wallcycle_stop(wcycle_, ewcWAIT_GPU_BONDED);
291
292     for (int fType : fTypesOnGpu)
293     {
294         if (fType != F_LJ14 && fType != F_COUL14)
295         {
296             enerd->term[fType] += vTot_[fType];
297         }
298     }
299
300     // Note: We do not support energy groups here
301     gmx_grppairener_t* grppener = &enerd->grpp;
302     GMX_RELEASE_ASSERT(grppener->nener == 1, "No energy group support for bondeds on the GPU");
303     grppener->ener[egLJ14][0] += vTot_[F_LJ14];
304     grppener->ener[egCOUL14][0] += vTot_[F_COUL14];
305 }
306
307 void GpuBonded::Impl::clearEnergies()
308 {
309     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
310     wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_GPU_BONDED);
311     clearDeviceBufferAsync(&d_vTot_, 0, F_NRE, deviceStream_);
312     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_BONDED);
313     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
314 }
315
316 // ---- GpuBonded
317
318 GpuBonded::GpuBonded(const gmx_ffparams_t& ffparams,
319                      const DeviceContext&  deviceContext,
320                      const DeviceStream&   deviceStream,
321                      gmx_wallcycle*        wcycle) :
322     impl_(new Impl(ffparams, deviceContext, deviceStream, wcycle))
323 {
324 }
325
326 GpuBonded::~GpuBonded() = default;
327
328 void GpuBonded::updateInteractionListsAndDeviceBuffers(ArrayRef<const int>           nbnxnAtomOrder,
329                                                        const InteractionDefinitions& idef,
330                                                        void*                         d_xq,
331                                                        DeviceBuffer<RVec>            d_f,
332                                                        DeviceBuffer<RVec>            d_fShift)
333 {
334     impl_->updateInteractionListsAndDeviceBuffers(nbnxnAtomOrder, idef, d_xq, d_f, d_fShift);
335 }
336
337 bool GpuBonded::haveInteractions() const
338 {
339     return impl_->haveInteractions();
340 }
341
342 void GpuBonded::launchEnergyTransfer()
343 {
344     impl_->launchEnergyTransfer();
345 }
346
347 void GpuBonded::waitAccumulateEnergyTerms(gmx_enerdata_t* enerd)
348 {
349     impl_->waitAccumulateEnergyTerms(enerd);
350 }
351
352 void GpuBonded::clearEnergies()
353 {
354     impl_->clearEnergies();
355 }
356
357 } // namespace gmx