dccfe1eed03ac6fa2ea38c4038921c90ca9f085b
[alexxy/gromacs.git] / src / gromacs / nbnxm / gpu_common.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,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  * \brief Common functions for the different NBNXN GPU implementations.
37  *
38  * \author Szilard Pall <pall.szilard@gmail.com>
39  *
40  * \ingroup module_nbnxm
41  */
42
43 #ifndef GMX_NBNXM_GPU_COMMON_H
44 #define GMX_NBNXM_GPU_COMMON_H
45
46 #include "config.h"
47
48 #include <string>
49
50 #if GMX_GPU_CUDA
51 #    include "cuda/nbnxm_cuda_types.h"
52 #endif
53
54 #if GMX_GPU_OPENCL
55 #    include "opencl/nbnxm_ocl_types.h"
56 #endif
57
58 #if GMX_GPU_SYCL
59 #    include "sycl/nbnxm_sycl_types.h"
60 #endif
61
62 #include "gromacs/gpu_utils/gpu_utils.h"
63 #include "gromacs/listed_forces/gpubonded.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/mdtypes/simulation_workload.h"
66 #include "gromacs/nbnxm/nbnxm.h"
67 #include "gromacs/pbcutil/ishift.h"
68 #include "gromacs/timing/gpu_timing.h"
69 #include "gromacs/timing/wallcycle.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/stringutil.h"
72
73 #include "gpu_common_utils.h"
74 #include "nbnxm_gpu.h"
75
76 namespace gmx
77 {
78 class GpuBonded;
79 }
80
81 namespace Nbnxm
82 {
83
84 /*! \brief Check that atom locality values are valid for the GPU module.
85  *
86  *  In the GPU module atom locality "all" is not supported, the local and
87  *  non-local ranges are treated separately.
88  *
89  *  \param[in] atomLocality atom locality specifier
90  */
91 static inline void validateGpuAtomLocality(const AtomLocality atomLocality)
92 {
93     std::string str = gmx::formatString(
94             "Invalid atom locality passed (%d); valid here is only "
95             "local (%d) or nonlocal (%d)",
96             static_cast<int>(atomLocality),
97             static_cast<int>(AtomLocality::Local),
98             static_cast<int>(AtomLocality::NonLocal));
99
100     GMX_ASSERT(atomLocality == AtomLocality::Local || atomLocality == AtomLocality::NonLocal, str.c_str());
101 }
102
103 /*! \brief Convert atom locality to interaction locality.
104  *
105  *  In the current implementation the this is straightforward conversion:
106  *  local to local, non-local to non-local.
107  *
108  *  \param[in] atomLocality Atom locality specifier
109  *  \returns                Interaction locality corresponding to the atom locality passed.
110  */
111 static inline InteractionLocality gpuAtomToInteractionLocality(const AtomLocality atomLocality)
112 {
113     validateGpuAtomLocality(atomLocality);
114
115     /* determine interaction locality from atom locality */
116     if (atomLocality == AtomLocality::Local)
117     {
118         return InteractionLocality::Local;
119     }
120     else if (atomLocality == AtomLocality::NonLocal)
121     {
122         return InteractionLocality::NonLocal;
123     }
124     else
125     {
126         gmx_incons("Wrong locality");
127     }
128 }
129
130
131 //NOLINTNEXTLINE(misc-definitions-in-headers)
132 void setupGpuShortRangeWork(NbnxmGpu* nb, const gmx::GpuBonded* gpuBonded, const gmx::InteractionLocality iLocality)
133 {
134     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
135
136     // There is short-range work if the pair list for the provided
137     // interaction locality contains entries or if there is any
138     // bonded work (as this is not split into local/nonlocal).
139     nb->haveWork[iLocality] = ((nb->plist[iLocality]->nsci != 0)
140                                || (gpuBonded != nullptr && gpuBonded->haveInteractions()));
141 }
142
143 /*! \brief Returns true if there is GPU short-range work for the given interaction locality.
144  *
145  * Note that as, unlike nonbonded tasks, bonded tasks are not split into local/nonlocal,
146  * and therefore if there are GPU offloaded bonded interactions, this function will return
147  * true for all interaction localities.
148  *
149  * \param[inout]  nb        Pointer to the nonbonded GPU data structure
150  * \param[in]     iLocality Interaction locality identifier
151  */
152 static bool haveGpuShortRangeWork(const NbnxmGpu& nb, const gmx::InteractionLocality iLocality)
153 {
154     return nb.haveWork[iLocality];
155 }
156
157 //NOLINTNEXTLINE(misc-definitions-in-headers)
158 bool haveGpuShortRangeWork(const NbnxmGpu* nb, const gmx::AtomLocality aLocality)
159 {
160     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
161
162     return haveGpuShortRangeWork(*nb, gpuAtomToInteractionLocality(aLocality));
163 }
164
165
166 /*! \brief Calculate atom range and return start index and length.
167  *
168  * \param[in] atomData Atom descriptor data structure
169  * \param[in] atomLocality Atom locality specifier
170  * \param[out] atomRangeBegin Starting index of the atom range in the atom data array.
171  * \param[out] atomRangeLen Atom range length in the atom data array.
172  */
173 template<typename AtomDataT>
174 static inline void getGpuAtomRange(const AtomDataT*   atomData,
175                                    const AtomLocality atomLocality,
176                                    int*               atomRangeBegin,
177                                    int*               atomRangeLen)
178 {
179     assert(atomData);
180     validateGpuAtomLocality(atomLocality);
181
182     /* calculate the atom data index range based on locality */
183     if (atomLocality == AtomLocality::Local)
184     {
185         *atomRangeBegin = 0;
186         *atomRangeLen   = atomData->numAtomsLocal;
187     }
188     else
189     {
190         *atomRangeBegin = atomData->numAtomsLocal;
191         *atomRangeLen   = atomData->numAtoms - atomData->numAtomsLocal;
192     }
193 }
194
195
196 /*! \brief Count pruning kernel time if either kernel has been triggered
197  *
198  *  We do the accounting for either of the two pruning kernel flavors:
199  *   - 1st pass prune: ran during the current step (prior to the force kernel);
200  *   - rolling prune:  ran at the end of the previous step (prior to the current step H2D xq);
201  *
202  * Note that the resetting of cu_timers_t::didPrune and cu_timers_t::didRollingPrune should happen
203  * after calling this function.
204  *
205  * \param[in] timers   structs with GPU timer objects
206  * \param[inout] timings  GPU task timing data
207  * \param[in] iloc        interaction locality
208  */
209 template<typename GpuTimers>
210 static void countPruneKernelTime(GpuTimers*                 timers,
211                                  gmx_wallclock_gpu_nbnxn_t* timings,
212                                  const InteractionLocality  iloc)
213 {
214     gpu_timers_t::Interaction& iTimers = timers->interaction[iloc];
215
216     // We might have not done any pruning (e.g. if we skipped with empty domains).
217     if (!iTimers.didPrune && !iTimers.didRollingPrune)
218     {
219         return;
220     }
221
222     if (iTimers.didPrune)
223     {
224         timings->pruneTime.c++;
225         timings->pruneTime.t += iTimers.prune_k.getLastRangeTime();
226     }
227
228     if (iTimers.didRollingPrune)
229     {
230         timings->dynamicPruneTime.c++;
231         timings->dynamicPruneTime.t += iTimers.rollingPrune_k.getLastRangeTime();
232     }
233 }
234
235 /*! \brief Reduce data staged internally in the nbnxn module.
236  *
237  * Shift forces and electrostatic/LJ energies copied from the GPU into
238  * a module-internal staging area are immediately reduced (CPU-side buffers passed)
239  * after having waited for the transfers' completion.
240  *
241  * Note that this function should always be called after the transfers into the
242  * staging buffers has completed.
243  *
244  * \param[in]  nbst           Nonbonded staging data
245  * \param[in]  iLocality      Interaction locality specifier
246  * \param[in]  reduceEnergies True if energy reduction should be done
247  * \param[in]  reduceFshift   True if shift force reduction should be done
248  * \param[out] e_lj           Variable to accumulate LJ energy into
249  * \param[out] e_el           Variable to accumulate electrostatic energy into
250  * \param[out] fshift         Pointer to the array of shift forces to accumulate into
251  */
252 static inline void gpu_reduce_staged_outputs(const NBStagingData&      nbst,
253                                              const InteractionLocality iLocality,
254                                              const bool                reduceEnergies,
255                                              const bool                reduceFshift,
256                                              real*                     e_lj,
257                                              real*                     e_el,
258                                              rvec*                     fshift)
259 {
260     /* add up energies and shift forces (only once at local F wait) */
261     if (iLocality == InteractionLocality::Local)
262     {
263         if (reduceEnergies)
264         {
265             *e_lj += *nbst.eLJ;
266             *e_el += *nbst.eElec;
267         }
268
269         if (reduceFshift)
270         {
271             for (int i = 0; i < SHIFTS; i++)
272             {
273                 rvec_inc(fshift[i], nbst.fShift[i]);
274             }
275         }
276     }
277 }
278
279 /*! \brief Do the per-step timing accounting of the nonbonded tasks.
280  *
281  *  Does timing accumulation and call-count increments for the nonbonded kernels.
282  *  Note that this function should be called after the current step's nonbonded
283  *  nonbonded tasks have completed with the exception of the rolling pruning kernels
284  *  that are accounted for during the following step.
285  *
286  * NOTE: if timing with multiple GPUs (streams) becomes possible, the
287  *      counters could end up being inconsistent due to not being incremented
288  *      on some of the node when this is skipped on empty local domains!
289  *
290  * \tparam     GpuTimers         GPU timers type
291  * \tparam     GpuPairlist       Pair list type
292  * \param[out] timings           Pointer to the NB GPU timings data
293  * \param[in]  timers            Pointer to GPU timers data
294  * \param[in]  plist             Pointer to the pair list data
295  * \param[in]  atomLocality      Atom locality specifier
296  * \param[in]  stepWork          Force schedule flags
297  * \param[in]  doTiming          True if timing is enabled.
298  *
299  */
300 template<typename GpuTimers, typename GpuPairlist>
301 static inline void gpu_accumulate_timings(gmx_wallclock_gpu_nbnxn_t* timings,
302                                           GpuTimers*                 timers,
303                                           const GpuPairlist*         plist,
304                                           AtomLocality               atomLocality,
305                                           const gmx::StepWorkload&   stepWork,
306                                           bool                       doTiming)
307 {
308     /* timing data accumulation */
309     if (!doTiming)
310     {
311         return;
312     }
313
314     /* determine interaction locality from atom locality */
315     const InteractionLocality iLocality        = gpuAtomToInteractionLocality(atomLocality);
316     const bool                didEnergyKernels = stepWork.computeEnergy;
317
318     /* only increase counter once (at local F wait) */
319     if (iLocality == InteractionLocality::Local)
320     {
321         timings->nb_c++;
322         timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].c += 1;
323     }
324
325     /* kernel timings */
326     timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].t +=
327             timers->interaction[iLocality].nb_k.getLastRangeTime();
328
329     /* X/q H2D and F D2H timings */
330     timings->nb_h2d_t += timers->xf[atomLocality].nb_h2d.getLastRangeTime();
331     timings->nb_d2h_t += timers->xf[atomLocality].nb_d2h.getLastRangeTime();
332
333     /* Count the pruning kernel times for both cases:1st pass (at search step)
334        and rolling pruning (if called at the previous step).
335        We do the accounting here as this is the only sync point where we
336        know (without checking or additional sync-ing) that prune tasks in
337        in the current stream have completed (having just blocking-waited
338        for the force D2H). */
339     countPruneKernelTime(timers, timings, iLocality);
340
341     /* only count atdat at pair-search steps (add only once, at local F wait) */
342     if (stepWork.doNeighborSearch && atomLocality == AtomLocality::Local)
343     {
344         /* atdat transfer timing */
345         timings->pl_h2d_c++;
346         timings->pl_h2d_t += timers->atdat.getLastRangeTime();
347     }
348
349     /* only count pair-list H2D when actually performed */
350     if (timers->interaction[iLocality].didPairlistH2D)
351     {
352         timings->pl_h2d_t += timers->interaction[iLocality].pl_h2d.getLastRangeTime();
353
354         /* Clear the timing flag for the next step */
355         timers->interaction[iLocality].didPairlistH2D = false;
356     }
357 }
358
359 /*! \brief Attempts to complete nonbonded GPU task.
360  *
361  * See documentation in nbnxm_gpu.h for details.
362  *
363  * \todo Move into shared source file, perhaps including
364  * cuda_runtime.h if needed for any remaining CUDA-specific
365  * objects.
366  */
367 //NOLINTNEXTLINE(misc-definitions-in-headers)
368 bool gpu_try_finish_task(NbnxmGpu*                nb,
369                          const gmx::StepWorkload& stepWork,
370                          const AtomLocality       aloc,
371                          real*                    e_lj,
372                          real*                    e_el,
373                          gmx::ArrayRef<gmx::RVec> shiftForces,
374                          GpuTaskCompletion        completionKind,
375                          gmx_wallcycle*           wcycle)
376 {
377     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
378
379     /* determine interaction locality from atom locality */
380     const InteractionLocality iLocality = gpuAtomToInteractionLocality(aloc);
381
382
383     // Transfers are launched and therefore need to be waited on if:
384     // - buffer ops is not offloaded
385     // - energies or virials are needed (on the local stream)
386     //
387     // (Note that useGpuFBufferOps and computeVirial are mutually exclusive
388     // in current code as virial steps do CPU reduction.)
389     const bool haveResultToWaitFor =
390             (!stepWork.useGpuFBufferOps
391              || (aloc == AtomLocality::Local && (stepWork.computeEnergy || stepWork.computeVirial)));
392
393     //  We skip when during the non-local phase there was actually no work to do.
394     //  This is consistent with nbnxn_gpu_launch_kernel but it also considers possible
395     //  bonded GPU work.
396     if ((iLocality == InteractionLocality::Local) || haveGpuShortRangeWork(*nb, iLocality))
397     {
398         // Query the state of the GPU stream and return early if we're not done
399         if (completionKind == GpuTaskCompletion::Check)
400         {
401             // To get the wcycle call count right, when in GpuTaskCompletion::Check mode,
402             // we start without counting and only when the task finished we issue a
403             // start/stop to increment.
404             // GpuTaskCompletion::Wait mode the timing is expected to be done in the caller.
405             wallcycle_start_nocount(wcycle, ewcWAIT_GPU_NB_L);
406
407             if (!haveStreamTasksCompleted(*nb->deviceStreams[iLocality]))
408             {
409                 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
410
411                 // Early return to skip the steps below that we have to do only
412                 // after the NB task completed
413                 return false;
414             }
415
416             wallcycle_increment_event_count(wcycle, ewcWAIT_GPU_NB_L);
417         }
418         else if (haveResultToWaitFor)
419         {
420             nb->deviceStreams[iLocality]->synchronize();
421         }
422
423         // TODO: this needs to be moved later because conditional wait could brake timing
424         // with a future OpenCL implementation, but with CUDA timing is anyway disabled
425         // in all cases where we skip the wait.
426         gpu_accumulate_timings(nb->timings, nb->timers, nb->plist[iLocality], aloc, stepWork, nb->bDoTime);
427
428         if (stepWork.computeEnergy || stepWork.computeVirial)
429         {
430             gpu_reduce_staged_outputs(nb->nbst,
431                                       iLocality,
432                                       stepWork.computeEnergy,
433                                       stepWork.computeVirial,
434                                       e_lj,
435                                       e_el,
436                                       as_rvec_array(shiftForces.data()));
437         }
438     }
439
440     /* Reset both pruning flags. */
441     if (nb->bDoTime)
442     {
443         nb->timers->interaction[iLocality].didPrune =
444                 nb->timers->interaction[iLocality].didRollingPrune = false;
445     }
446
447     /* Turn off initial list pruning (doesn't hurt if this is not pair-search step). */
448     nb->plist[iLocality]->haveFreshList = false;
449
450     return true;
451 }
452
453 /*! \brief
454  * Wait for the asynchronously launched nonbonded tasks and data
455  * transfers to finish.
456  *
457  * Also does timing accounting and reduction of the internal staging buffers.
458  * As this is called at the end of the step, it also resets the pair list and
459  * pruning flags.
460  *
461  * \param[in] nb The nonbonded data GPU structure
462  * \param[in]  stepWork     Force schedule flags
463  * \param[in] aloc Atom locality identifier
464  * \param[out] e_lj Pointer to the LJ energy output to accumulate into
465  * \param[out] e_el Pointer to the electrostatics energy output to accumulate into
466  * \param[out] shiftForces Shift forces buffer to accumulate into
467  * \param[out] wcycle Pointer to wallcycle data structure
468  * \return            The number of cycles the gpu wait took
469  */
470 //NOLINTNEXTLINE(misc-definitions-in-headers) TODO: move into source file
471 float gpu_wait_finish_task(NbnxmGpu*                nb,
472                            const gmx::StepWorkload& stepWork,
473                            AtomLocality             aloc,
474                            real*                    e_lj,
475                            real*                    e_el,
476                            gmx::ArrayRef<gmx::RVec> shiftForces,
477                            gmx_wallcycle*           wcycle)
478 {
479     auto cycleCounter = (gpuAtomToInteractionLocality(aloc) == InteractionLocality::Local)
480                                 ? ewcWAIT_GPU_NB_L
481                                 : ewcWAIT_GPU_NB_NL;
482
483     wallcycle_start(wcycle, cycleCounter);
484     gpu_try_finish_task(nb, stepWork, aloc, e_lj, e_el, shiftForces, GpuTaskCompletion::Wait, wcycle);
485     float waitTime = wallcycle_stop(wcycle, cycleCounter);
486
487     return waitTime;
488 }
489
490 } // namespace Nbnxm
491
492 #endif