2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2017,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.
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.
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.
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.
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.
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.
36 * \brief Common functions for the different NBNXN GPU implementations.
38 * \author Szilard Pall <pall.szilard@gmail.com>
40 * \ingroup module_nbnxm
43 #ifndef GMX_NBNXM_GPU_COMMON_H
44 #define GMX_NBNXM_GPU_COMMON_H
50 #if GMX_GPU == GMX_GPU_CUDA
51 #include "cuda/nbnxm_cuda_types.h"
54 #if GMX_GPU == GMX_GPU_OPENCL
55 #include "opencl/nbnxm_ocl_types.h"
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"
69 #include "gpu_common_utils.h"
70 #include "nbnxm_gpu.h"
80 /*! \brief Check that atom locality values are valid for the GPU module.
82 * In the GPU module atom locality "all" is not supported, the local and
83 * non-local ranges are treated separately.
85 * \param[in] atomLocality atom locality specifier
88 validateGpuAtomLocality(const AtomLocality atomLocality)
90 std::string str = gmx::formatString("Invalid atom locality passed (%d); valid here is only "
91 "local (%d) or nonlocal (%d)",
92 static_cast<int>(atomLocality),
93 static_cast<int>(AtomLocality::Local),
94 static_cast<int>(AtomLocality::NonLocal));
96 GMX_ASSERT(atomLocality == AtomLocality::Local || atomLocality == AtomLocality::NonLocal, str.c_str());
99 /*! \brief Convert atom locality to interaction locality.
101 * In the current implementation the this is straightforward conversion:
102 * local to local, non-local to non-local.
104 * \param[in] atomLocality Atom locality specifier
105 * \returns Interaction locality corresponding to the atom locality passed.
107 static inline InteractionLocality
108 gpuAtomToInteractionLocality(const AtomLocality atomLocality)
110 validateGpuAtomLocality(atomLocality);
112 /* determine interaction locality from atom locality */
113 if (atomLocality == AtomLocality::Local)
115 return InteractionLocality::Local;
117 else if (atomLocality == AtomLocality::NonLocal)
119 return InteractionLocality::NonLocal;
123 gmx_incons("Wrong locality");
128 //NOLINTNEXTLINE(misc-definitions-in-headers)
129 void setupGpuShortRangeWork(gmx_nbnxn_gpu_t *nb,
130 const gmx::GpuBonded *gpuBonded,
131 const gmx::InteractionLocality iLocality)
133 GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
135 // There is short-range work if the pair list for the provided
136 // interaction locality contains entries or if there is any
137 // bonded work (as this is not split into local/nonlocal).
138 nb->haveWork[iLocality] =
139 ((nb->plist[iLocality]->nsci != 0) ||
140 (gpuBonded != nullptr && gpuBonded->haveInteractions()));
143 /*! \brief Returns true if there is GPU short-range work for the given interaction locality.
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.
149 * \param[inout] nb Pointer to the nonbonded GPU data structure
150 * \param[in] iLocality Interaction locality identifier
153 haveGpuShortRangeWork(const gmx_nbnxn_gpu_t &nb,
154 const gmx::InteractionLocality iLocality)
156 return nb.haveWork[iLocality];
159 //NOLINTNEXTLINE(misc-definitions-in-headers)
160 bool haveGpuShortRangeWork(const gmx_nbnxn_gpu_t *nb,
161 const gmx::AtomLocality aLocality)
163 GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
165 return haveGpuShortRangeWork(*nb, gpuAtomToInteractionLocality(aLocality));
170 /*! \brief Calculate atom range and return start index and length.
172 * \param[in] atomData Atom descriptor data structure
173 * \param[in] atomLocality Atom locality specifier
174 * \param[out] atomRangeBegin Starting index of the atom range in the atom data array.
175 * \param[out] atomRangeLen Atom range length in the atom data array.
177 template <typename AtomDataT>
179 getGpuAtomRange(const AtomDataT *atomData,
180 const AtomLocality atomLocality,
185 validateGpuAtomLocality(atomLocality);
187 /* calculate the atom data index range based on locality */
188 if (atomLocality == AtomLocality::Local)
191 *atomRangeLen = atomData->natoms_local;
195 *atomRangeBegin = atomData->natoms_local;
196 *atomRangeLen = atomData->natoms - atomData->natoms_local;
201 /*! \brief Count pruning kernel time if either kernel has been triggered
203 * We do the accounting for either of the two pruning kernel flavors:
204 * - 1st pass prune: ran during the current step (prior to the force kernel);
205 * - rolling prune: ran at the end of the previous step (prior to the current step H2D xq);
207 * Note that the resetting of cu_timers_t::didPrune and cu_timers_t::didRollingPrune should happen
208 * after calling this function.
210 * \param[in] timers structs with GPU timer objects
211 * \param[inout] timings GPU task timing data
212 * \param[in] iloc interaction locality
214 template <typename GpuTimers>
215 static void countPruneKernelTime(GpuTimers *timers,
216 gmx_wallclock_gpu_nbnxn_t *timings,
217 const InteractionLocality iloc)
219 gpu_timers_t::Interaction &iTimers = timers->interaction[iloc];
221 // We might have not done any pruning (e.g. if we skipped with empty domains).
222 if (!iTimers.didPrune &&
223 !iTimers.didRollingPrune)
228 if (iTimers.didPrune)
230 timings->pruneTime.c++;
231 timings->pruneTime.t += iTimers.prune_k.getLastRangeTime();
234 if (iTimers.didRollingPrune)
236 timings->dynamicPruneTime.c++;
237 timings->dynamicPruneTime.t += iTimers.rollingPrune_k.getLastRangeTime();
241 /*! \brief Reduce data staged internally in the nbnxn module.
243 * Shift forces and electrostatic/LJ energies copied from the GPU into
244 * a module-internal staging area are immediately reduced (CPU-side buffers passed)
245 * after having waited for the transfers' completion.
247 * Note that this function should always be called after the transfers into the
248 * staging buffers has completed.
250 * \tparam StagingData Type of staging data
251 * \param[in] nbst Nonbonded staging data
252 * \param[in] iLocality Interaction locality specifier
253 * \param[in] reduceEnergies True if energy reduction should be done
254 * \param[in] reduceFshift True if shift force reduction should be done
255 * \param[out] e_lj Variable to accumulate LJ energy into
256 * \param[out] e_el Variable to accumulate electrostatic energy into
257 * \param[out] fshift Pointer to the array of shift forces to accumulate into
259 template <typename StagingData>
261 gpu_reduce_staged_outputs(const StagingData &nbst,
262 const InteractionLocality iLocality,
263 const bool reduceEnergies,
264 const bool reduceFshift,
269 /* add up energies and shift forces (only once at local F wait) */
270 if (iLocality == InteractionLocality::Local)
280 for (int i = 0; i < SHIFTS; i++)
282 rvec_inc(fshift[i], nbst.fshift[i]);
288 /*! \brief Do the per-step timing accounting of the nonbonded tasks.
290 * Does timing accumulation and call-count increments for the nonbonded kernels.
291 * Note that this function should be called after the current step's nonbonded
292 * nonbonded tasks have completed with the exception of the rolling pruning kernels
293 * that are accounted for during the following step.
295 * NOTE: if timing with multiple GPUs (streams) becomes possible, the
296 * counters could end up being inconsistent due to not being incremented
297 * on some of the node when this is skipped on empty local domains!
299 * \tparam GpuTimers GPU timers type
300 * \tparam GpuPairlist Pair list type
301 * \param[out] timings Pointer to the NB GPU timings data
302 * \param[in] timers Pointer to GPU timers data
303 * \param[in] plist Pointer to the pair list data
304 * \param[in] atomLocality Atom locality specifier
305 * \param[in] stepWork Force schedule flags
306 * \param[in] doTiming True if timing is enabled.
309 template <typename GpuTimers, typename GpuPairlist>
311 gpu_accumulate_timings(gmx_wallclock_gpu_nbnxn_t *timings,
313 const GpuPairlist *plist,
314 AtomLocality atomLocality,
315 const gmx::StepWorkload &stepWork,
318 /* timing data accumulation */
324 /* determine interaction locality from atom locality */
325 const InteractionLocality iLocality = gpuAtomToInteractionLocality(atomLocality);
326 const bool didEnergyKernels = stepWork.computeEnergy;
328 /* only increase counter once (at local F wait) */
329 if (iLocality == InteractionLocality::Local)
332 timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].c += 1;
336 timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].t +=
337 timers->interaction[iLocality].nb_k.getLastRangeTime();
339 /* X/q H2D and F D2H timings */
340 timings->nb_h2d_t += timers->xf[atomLocality].nb_h2d.getLastRangeTime();
341 timings->nb_d2h_t += timers->xf[atomLocality].nb_d2h.getLastRangeTime();
343 /* Count the pruning kernel times for both cases:1st pass (at search step)
344 and rolling pruning (if called at the previous step).
345 We do the accounting here as this is the only sync point where we
346 know (without checking or additional sync-ing) that prune tasks in
347 in the current stream have completed (having just blocking-waited
348 for the force D2H). */
349 countPruneKernelTime(timers, timings, iLocality);
351 /* only count atdat at pair-search steps (add only once, at local F wait) */
352 if (stepWork.doNeighborSearch && atomLocality == AtomLocality::Local)
354 /* atdat transfer timing */
356 timings->pl_h2d_t += timers->atdat.getLastRangeTime();
359 /* only count pair-list H2D when actually performed */
360 if (timers->interaction[iLocality].didPairlistH2D)
362 timings->pl_h2d_t += timers->interaction[iLocality].pl_h2d.getLastRangeTime();
364 /* Clear the timing flag for the next step */
365 timers->interaction[iLocality].didPairlistH2D = false;
369 /*! \brief Attempts to complete nonbonded GPU task.
371 * See documentation in nbnxm_gpu.h for details.
373 * \todo Move into shared source file with gmx_compile_cpp_as_cuda
375 //NOLINTNEXTLINE(misc-definitions-in-headers)
376 bool gpu_try_finish_task(gmx_nbnxn_gpu_t *nb,
377 const gmx::StepWorkload &stepWork,
378 const AtomLocality aloc,
381 gmx::ArrayRef<gmx::RVec> shiftForces,
382 GpuTaskCompletion completionKind,
383 gmx_wallcycle *wcycle)
385 GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
387 /* determine interaction locality from atom locality */
388 const InteractionLocality iLocality = gpuAtomToInteractionLocality(aloc);
391 // Transfers are launched and therefore need to be waited on if:
392 // - buffer ops is not offloaded
393 // - energies or virials are needed (on the local stream)
395 // (Note that useGpuFBufferOps and computeVirial are mutually exclusive
396 // in current code as virial steps do CPU reduction.)
397 const bool haveResultToWaitFor =
398 (!stepWork.useGpuFBufferOps ||
399 (aloc == AtomLocality::Local && (stepWork.computeEnergy || stepWork.computeVirial)));
401 // We skip when during the non-local phase there was actually no work to do.
402 // This is consistent with nbnxn_gpu_launch_kernel but it also considers possible
404 if ((iLocality == InteractionLocality::Local) || haveGpuShortRangeWork(*nb, iLocality))
406 // Query the state of the GPU stream and return early if we're not done
407 if (completionKind == GpuTaskCompletion::Check)
409 // To get the wcycle call count right, when in GpuTaskCompletion::Check mode,
410 // we start without counting and only when the task finished we issue a
411 // start/stop to increment.
412 // GpuTaskCompletion::Wait mode the timing is expected to be done in the caller.
413 wallcycle_start_nocount(wcycle, ewcWAIT_GPU_NB_L);
415 if (!haveStreamTasksCompleted(nb->stream[iLocality]))
417 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
419 // Early return to skip the steps below that we have to do only
420 // after the NB task completed
424 wallcycle_increment_event_count(wcycle, ewcWAIT_GPU_NB_L);
426 else if (haveResultToWaitFor)
428 gpuStreamSynchronize(nb->stream[iLocality]);
431 // TODO: this needs to be moved later because conditional wait could brake timing
432 // with a future OpenCL implementation, but with CUDA timing is anyway disabled
433 // in all cases where we skip the wait.
434 gpu_accumulate_timings(nb->timings, nb->timers, nb->plist[iLocality], aloc, stepWork,
437 if (stepWork.computeEnergy || stepWork.computeVirial)
439 gpu_reduce_staged_outputs(nb->nbst, iLocality, stepWork.computeEnergy, stepWork.computeVirial,
440 e_lj, e_el, as_rvec_array(shiftForces.data()));
444 /* Always reset both pruning flags (doesn't hurt doing it even when timing is off). */
445 nb->timers->interaction[iLocality].didPrune = nb->timers->interaction[iLocality].didRollingPrune = false;
447 /* Turn off initial list pruning (doesn't hurt if this is not pair-search step). */
448 nb->plist[iLocality]->haveFreshList = false;
454 * Wait for the asynchronously launched nonbonded tasks and data
455 * transfers to finish.
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
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
470 //NOLINTNEXTLINE(misc-definitions-in-headers) TODO: move into source file
471 float gpu_wait_finish_task(gmx_nbnxn_gpu_t *nb,
472 const gmx::StepWorkload &stepWork,
476 gmx::ArrayRef<gmx::RVec> shiftForces,
477 gmx_wallcycle *wcycle)
480 (gpuAtomToInteractionLocality(aloc) == InteractionLocality::Local) ? ewcWAIT_GPU_NB_L : ewcWAIT_GPU_NB_NL;
482 wallcycle_start(wcycle, cycleCounter);
483 gpu_try_finish_task(nb, stepWork, aloc, e_lj, e_el, shiftForces,
484 GpuTaskCompletion::Wait, wcycle);
485 float waitTime = wallcycle_stop(wcycle, cycleCounter);