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