Rename GpuBonded into ListedForcesGpu
[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/math/vec.h"
64 #include "gromacs/mdtypes/simulation_workload.h"
65 #include "gromacs/nbnxm/nbnxm.h"
66 #include "gromacs/pbcutil/ishift.h"
67 #include "gromacs/timing/gpu_timing.h"
68 #include "gromacs/timing/wallcycle.h"
69 #include "gromacs/utility/stringutil.h"
70
71 #include "gpu_common_utils.h"
72 #include "nbnxm_gpu.h"
73
74 namespace gmx
75 {
76 class ListedForcesGpu;
77 }
78
79 namespace Nbnxm
80 {
81
82 /*! \brief Count pruning kernel time if either kernel has been triggered
83  *
84  *  We do the accounting for either of the two pruning kernel flavors:
85  *   - 1st pass prune: ran during the current step (prior to the force kernel);
86  *   - rolling prune:  ran at the end of the previous step (prior to the current step H2D xq);
87  *
88  * Note that the resetting of Nbnxm::GpuTimers::didPrune and Nbnxm::GpuTimers::didRollingPrune
89  * should happen after calling this function.
90  *
91  * \param[in] timers   structs with GPU timer objects
92  * \param[inout] timings  GPU task timing data
93  * \param[in] iloc        interaction locality
94  */
95 static void countPruneKernelTime(Nbnxm::GpuTimers*          timers,
96                                  gmx_wallclock_gpu_nbnxn_t* timings,
97                                  const InteractionLocality  iloc)
98 {
99     GpuTimers::Interaction& iTimers = timers->interaction[iloc];
100
101     // We might have not done any pruning (e.g. if we skipped with empty domains).
102     if (!iTimers.didPrune && !iTimers.didRollingPrune)
103     {
104         return;
105     }
106
107     if (iTimers.didPrune)
108     {
109         timings->pruneTime.c++;
110         timings->pruneTime.t += iTimers.prune_k.getLastRangeTime();
111     }
112
113     if (iTimers.didRollingPrune)
114     {
115         timings->dynamicPruneTime.c++;
116         timings->dynamicPruneTime.t += iTimers.rollingPrune_k.getLastRangeTime();
117     }
118 }
119
120 /*! \brief Reduce data staged internally in the nbnxn module.
121  *
122  * Shift forces and electrostatic/LJ energies copied from the GPU into
123  * a module-internal staging area are immediately reduced (CPU-side buffers passed)
124  * after having waited for the transfers' completion.
125  *
126  * Note that this function should always be called after the transfers into the
127  * staging buffers has completed.
128  *
129  * \param[in]  nbst           Nonbonded staging data
130  * \param[in]  iLocality      Interaction locality specifier
131  * \param[in]  reduceEnergies True if energy reduction should be done
132  * \param[in]  reduceFshift   True if shift force reduction should be done
133  * \param[out] e_lj           Variable to accumulate LJ energy into
134  * \param[out] e_el           Variable to accumulate electrostatic energy into
135  * \param[out] fshift         Pointer to the array of shift forces to accumulate into
136  */
137 static inline void gpu_reduce_staged_outputs(const NBStagingData&      nbst,
138                                              const InteractionLocality iLocality,
139                                              const bool                reduceEnergies,
140                                              const bool                reduceFshift,
141                                              real*                     e_lj,
142                                              real*                     e_el,
143                                              rvec*                     fshift)
144 {
145     /* add up energies and shift forces (only once at local F wait) */
146     if (iLocality == InteractionLocality::Local)
147     {
148         if (reduceEnergies)
149         {
150             *e_lj += *nbst.eLJ;
151             *e_el += *nbst.eElec;
152         }
153
154         if (reduceFshift)
155         {
156             for (int i = 0; i < gmx::c_numShiftVectors; i++)
157             {
158                 rvec_inc(fshift[i], nbst.fShift[i]);
159             }
160         }
161     }
162 }
163
164 /*! \brief Do the per-step timing accounting of the nonbonded tasks.
165  *
166  *  Does timing accumulation and call-count increments for the nonbonded kernels.
167  *  Note that this function should be called after the current step's nonbonded
168  *  nonbonded tasks have completed with the exception of the rolling pruning kernels
169  *  that are accounted for during the following step.
170  *
171  * NOTE: if timing with multiple GPUs (streams) becomes possible, the
172  *      counters could end up being inconsistent due to not being incremented
173  *      on some of the node when this is skipped on empty local domains!
174  *
175  * \tparam     GpuPairlist       Pair list type
176  * \param[out] timings           Pointer to the NB GPU timings data
177  * \param[in]  timers            Pointer to GPU timers data
178  * \param[in]  plist             Pointer to the pair list data
179  * \param[in]  atomLocality      Atom locality specifier
180  * \param[in]  stepWork          Force schedule flags
181  * \param[in]  doTiming          True if timing is enabled.
182  *
183  */
184 template<typename GpuPairlist>
185 static inline void gpu_accumulate_timings(gmx_wallclock_gpu_nbnxn_t* timings,
186                                           Nbnxm::GpuTimers*          timers,
187                                           const GpuPairlist*         plist,
188                                           AtomLocality               atomLocality,
189                                           const gmx::StepWorkload&   stepWork,
190                                           bool                       doTiming)
191 {
192     /* timing data accumulation */
193     if (!doTiming)
194     {
195         return;
196     }
197
198     /* determine interaction locality from atom locality */
199     const InteractionLocality iLocality        = atomToInteractionLocality(atomLocality);
200     const bool                didEnergyKernels = stepWork.computeEnergy;
201
202     /* only increase counter once (at local F wait) */
203     if (iLocality == InteractionLocality::Local)
204     {
205         timings->nb_c++;
206         timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].c += 1;
207     }
208
209     /* kernel timings */
210     timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].t +=
211             timers->interaction[iLocality].nb_k.getLastRangeTime();
212
213     /* X/q H2D and F D2H timings */
214     timings->nb_h2d_t += timers->xf[atomLocality].nb_h2d.getLastRangeTime();
215     timings->nb_d2h_t += timers->xf[atomLocality].nb_d2h.getLastRangeTime();
216
217     /* Count the pruning kernel times for both cases:1st pass (at search step)
218        and rolling pruning (if called at the previous step).
219        We do the accounting here as this is the only sync point where we
220        know (without checking or additional sync-ing) that prune tasks in
221        in the current stream have completed (having just blocking-waited
222        for the force D2H). */
223     countPruneKernelTime(timers, timings, iLocality);
224
225     /* only count atdat at pair-search steps (add only once, at local F wait) */
226     if (stepWork.doNeighborSearch && atomLocality == AtomLocality::Local)
227     {
228         /* atdat transfer timing */
229         timings->pl_h2d_c++;
230         timings->pl_h2d_t += timers->atdat.getLastRangeTime();
231     }
232
233     /* only count pair-list H2D when actually performed */
234     if (timers->interaction[iLocality].didPairlistH2D)
235     {
236         timings->pl_h2d_t += timers->interaction[iLocality].pl_h2d.getLastRangeTime();
237
238         /* Clear the timing flag for the next step */
239         timers->interaction[iLocality].didPairlistH2D = false;
240     }
241 }
242
243 /*! \brief Attempts to complete nonbonded GPU task.
244  *
245  * See documentation in nbnxm_gpu.h for details.
246  *
247  * \todo Move into shared source file, perhaps including
248  * cuda_runtime.h if needed for any remaining CUDA-specific
249  * objects.
250  */
251 //NOLINTNEXTLINE(misc-definitions-in-headers)
252 bool gpu_try_finish_task(NbnxmGpu*                nb,
253                          const gmx::StepWorkload& stepWork,
254                          const AtomLocality       aloc,
255                          real*                    e_lj,
256                          real*                    e_el,
257                          gmx::ArrayRef<gmx::RVec> shiftForces,
258                          GpuTaskCompletion        completionKind,
259                          gmx_wallcycle*           wcycle)
260 {
261     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
262
263     /* determine interaction locality from atom locality */
264     const InteractionLocality iLocality = atomToInteractionLocality(aloc);
265
266
267     // Transfers are launched and therefore need to be waited on if:
268     // - buffer ops is not offloaded
269     // - energies or virials are needed (on the local stream)
270     //
271     // (Note that useGpuFBufferOps and computeVirial are mutually exclusive
272     // in current code as virial steps do CPU reduction.)
273     const bool haveResultToWaitFor =
274             (!stepWork.useGpuFBufferOps
275              || (aloc == AtomLocality::Local && (stepWork.computeEnergy || stepWork.computeVirial)));
276
277     //  We skip when during the non-local phase there was actually no work to do.
278     //  This is consistent with nbnxn_gpu_launch_kernel but it also considers possible
279     //  bonded GPU work.
280     if ((iLocality == InteractionLocality::Local) || haveGpuShortRangeWork(nb, iLocality))
281     {
282         // Query the state of the GPU stream and return early if we're not done
283         if (completionKind == GpuTaskCompletion::Check)
284         {
285             // To get the wcycle call count right, when in GpuTaskCompletion::Check mode,
286             // we start without counting and only when the task finished we issue a
287             // start/stop to increment.
288             // GpuTaskCompletion::Wait mode the timing is expected to be done in the caller.
289             wallcycle_start_nocount(wcycle, WallCycleCounter::WaitGpuNbL);
290
291             if (!haveStreamTasksCompleted(*nb->deviceStreams[iLocality]))
292             {
293                 wallcycle_stop(wcycle, WallCycleCounter::WaitGpuNbL);
294
295                 // Early return to skip the steps below that we have to do only
296                 // after the NB task completed
297                 return false;
298             }
299
300             wallcycle_increment_event_count(wcycle, WallCycleCounter::WaitGpuNbL);
301         }
302         else if (haveResultToWaitFor)
303         {
304             nb->deviceStreams[iLocality]->synchronize();
305         }
306
307         // TODO: this needs to be moved later because conditional wait could brake timing
308         // with a future OpenCL implementation, but with CUDA timing is anyway disabled
309         // in all cases where we skip the wait.
310         gpu_accumulate_timings(nb->timings, nb->timers, nb->plist[iLocality], aloc, stepWork, nb->bDoTime);
311
312         if (stepWork.computeEnergy || stepWork.computeVirial)
313         {
314             gpu_reduce_staged_outputs(nb->nbst,
315                                       iLocality,
316                                       stepWork.computeEnergy,
317                                       stepWork.computeVirial,
318                                       e_lj,
319                                       e_el,
320                                       as_rvec_array(shiftForces.data()));
321         }
322     }
323
324     /* Reset both pruning flags. */
325     if (nb->bDoTime)
326     {
327         nb->timers->interaction[iLocality].didPrune =
328                 nb->timers->interaction[iLocality].didRollingPrune = false;
329     }
330
331     /* Turn off initial list pruning (doesn't hurt if this is not pair-search step). */
332     nb->plist[iLocality]->haveFreshList = false;
333
334     return true;
335 }
336
337 /*! \brief
338  * Wait for the asynchronously launched nonbonded tasks and data
339  * transfers to finish.
340  *
341  * Also does timing accounting and reduction of the internal staging buffers.
342  * As this is called at the end of the step, it also resets the pair list and
343  * pruning flags.
344  *
345  * \param[in] nb The nonbonded data GPU structure
346  * \param[in]  stepWork     Force schedule flags
347  * \param[in] aloc Atom locality identifier
348  * \param[out] e_lj Pointer to the LJ energy output to accumulate into
349  * \param[out] e_el Pointer to the electrostatics energy output to accumulate into
350  * \param[out] shiftForces Shift forces buffer to accumulate into
351  * \param[out] wcycle Pointer to wallcycle data structure
352  * \return            The number of cycles the gpu wait took
353  */
354 //NOLINTNEXTLINE(misc-definitions-in-headers) TODO: move into source file
355 float gpu_wait_finish_task(NbnxmGpu*                nb,
356                            const gmx::StepWorkload& stepWork,
357                            AtomLocality             aloc,
358                            real*                    e_lj,
359                            real*                    e_el,
360                            gmx::ArrayRef<gmx::RVec> shiftForces,
361                            gmx_wallcycle*           wcycle)
362 {
363     auto cycleCounter = (atomToInteractionLocality(aloc) == InteractionLocality::Local)
364                                 ? WallCycleCounter::WaitGpuNbL
365                                 : WallCycleCounter::WaitGpuNbNL;
366
367     wallcycle_start(wcycle, cycleCounter);
368     gpu_try_finish_task(nb, stepWork, aloc, e_lj, e_el, shiftForces, GpuTaskCompletion::Wait, wcycle);
369     float waitTime = wallcycle_stop(wcycle, cycleCounter);
370
371     return waitTime;
372 }
373
374 } // namespace Nbnxm
375
376 #endif