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