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