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