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/math/vec.h"
60 #include "gromacs/mdlib/force_flags.h"
61 #include "gromacs/nbnxm/nbnxm.h"
62 #include "gromacs/pbcutil/ishift.h"
63 #include "gromacs/timing/gpu_timing.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/stringutil.h"
67 #include "gpu_common_utils.h"
68 #include "nbnxm_gpu.h"
73 /*! \brief Check that atom locality values are valid for the GPU module.
75 * In the GPU module atom locality "all" is not supported, the local and
76 * non-local ranges are treated separately.
78 * \param[in] atomLocality atom locality specifier
81 validateGpuAtomLocality(const AtomLocality atomLocality)
83 std::string str = gmx::formatString("Invalid atom locality passed (%d); valid here is only "
84 "local (%d) or nonlocal (%d)",
85 static_cast<int>(atomLocality),
86 static_cast<int>(AtomLocality::Local),
87 static_cast<int>(AtomLocality::NonLocal));
89 GMX_ASSERT(atomLocality == AtomLocality::Local || atomLocality == AtomLocality::NonLocal, str.c_str());
92 /*! \brief Convert atom locality to interaction locality.
94 * In the current implementation the this is straightforward conversion:
95 * local to local, non-local to non-local.
97 * \param[in] atomLocality Atom locality specifier
98 * \returns Interaction locality corresponding to the atom locality passed.
100 static inline InteractionLocality
101 gpuAtomToInteractionLocality(const AtomLocality atomLocality)
103 validateGpuAtomLocality(atomLocality);
105 /* determine interaction locality from atom locality */
106 if (atomLocality == AtomLocality::Local)
108 return InteractionLocality::Local;
110 else if (atomLocality == AtomLocality::NonLocal)
112 return InteractionLocality::NonLocal;
116 gmx_incons("Wrong locality");
120 /*! \brief Calculate atom range and return start index and length.
122 * \param[in] atomData Atom descriptor data structure
123 * \param[in] atomLocality Atom locality specifier
124 * \param[out] atomRangeBegin Starting index of the atom range in the atom data array.
125 * \param[out] atomRangeLen Atom range length in the atom data array.
127 template <typename AtomDataT>
129 getGpuAtomRange(const AtomDataT *atomData,
130 const AtomLocality atomLocality,
135 validateGpuAtomLocality(atomLocality);
137 /* calculate the atom data index range based on locality */
138 if (atomLocality == AtomLocality::Local)
141 *atomRangeLen = atomData->natoms_local;
145 *atomRangeBegin = atomData->natoms_local;
146 *atomRangeLen = atomData->natoms - atomData->natoms_local;
151 /*! \brief Count pruning kernel time if either kernel has been triggered
153 * We do the accounting for either of the two pruning kernel flavors:
154 * - 1st pass prune: ran during the current step (prior to the force kernel);
155 * - rolling prune: ran at the end of the previous step (prior to the current step H2D xq);
157 * Note that the resetting of cu_timers_t::didPrune and cu_timers_t::didRollingPrune should happen
158 * after calling this function.
160 * \param[in] timers structs with GPU timer objects
161 * \param[inout] timings GPU task timing data
162 * \param[in] iloc interaction locality
164 template <typename GpuTimers>
165 static void countPruneKernelTime(GpuTimers *timers,
166 gmx_wallclock_gpu_nbnxn_t *timings,
167 const InteractionLocality iloc)
169 gpu_timers_t::Interaction &iTimers = timers->interaction[iloc];
171 // We might have not done any pruning (e.g. if we skipped with empty domains).
172 if (!iTimers.didPrune &&
173 !iTimers.didRollingPrune)
178 if (iTimers.didPrune)
180 timings->pruneTime.c++;
181 timings->pruneTime.t += iTimers.prune_k.getLastRangeTime();
184 if (iTimers.didRollingPrune)
186 timings->dynamicPruneTime.c++;
187 timings->dynamicPruneTime.t += iTimers.rollingPrune_k.getLastRangeTime();
191 /*! \brief Reduce data staged internally in the nbnxn module.
193 * Shift forces and electrostatic/LJ energies copied from the GPU into
194 * a module-internal staging area are immediately reduced (CPU-side buffers passed)
195 * after having waited for the transfers' completion.
197 * Note that this function should always be called after the transfers into the
198 * staging buffers has completed.
200 * \tparam StagingData Type of staging data
201 * \param[in] nbst Nonbonded staging data
202 * \param[in] iLocality Interaction locality specifier
203 * \param[in] reduceEnergies True if energy reduction should be done
204 * \param[in] reduceFshift True if shift force reduction should be done
205 * \param[out] e_lj Variable to accumulate LJ energy into
206 * \param[out] e_el Variable to accumulate electrostatic energy into
207 * \param[out] fshift Pointer to the array of shift forces to accumulate into
209 template <typename StagingData>
211 gpu_reduce_staged_outputs(const StagingData &nbst,
212 const InteractionLocality iLocality,
213 const bool reduceEnergies,
214 const bool reduceFshift,
219 /* add up energies and shift forces (only once at local F wait) */
220 if (iLocality == InteractionLocality::Local)
230 for (int i = 0; i < SHIFTS; i++)
232 rvec_inc(fshift[i], nbst.fshift[i]);
238 /*! \brief Do the per-step timing accounting of the nonbonded tasks.
240 * Does timing accumulation and call-count increments for the nonbonded kernels.
241 * Note that this function should be called after the current step's nonbonded
242 * nonbonded tasks have completed with the exception of the rolling pruning kernels
243 * that are accounted for during the following step.
245 * NOTE: if timing with multiple GPUs (streams) becomes possible, the
246 * counters could end up being inconsistent due to not being incremented
247 * on some of the node when this is skipped on empty local domains!
249 * \tparam GpuTimers GPU timers type
250 * \tparam GpuPairlist Pair list type
251 * \param[out] timings Pointer to the NB GPU timings data
252 * \param[in] timers Pointer to GPU timers data
253 * \param[in] plist Pointer to the pair list data
254 * \param[in] atomLocality Atom locality specifier
255 * \param[in] didEnergyKernels True if energy kernels have been called in the current step
256 * \param[in] doTiming True if timing is enabled.
259 template <typename GpuTimers, typename GpuPairlist>
261 gpu_accumulate_timings(gmx_wallclock_gpu_nbnxn_t *timings,
263 const GpuPairlist *plist,
264 AtomLocality atomLocality,
265 bool didEnergyKernels,
268 /* timing data accumulation */
274 /* determine interaction locality from atom locality */
275 const InteractionLocality iLocality = gpuAtomToInteractionLocality(atomLocality);
277 /* only increase counter once (at local F wait) */
278 if (iLocality == InteractionLocality::Local)
281 timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].c += 1;
285 timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].t +=
286 timers->interaction[iLocality].nb_k.getLastRangeTime();
288 /* X/q H2D and F D2H timings */
289 timings->nb_h2d_t += timers->xf[atomLocality].nb_h2d.getLastRangeTime();
290 timings->nb_d2h_t += timers->xf[atomLocality].nb_d2h.getLastRangeTime();
292 /* Count the pruning kernel times for both cases:1st pass (at search step)
293 and rolling pruning (if called at the previous step).
294 We do the accounting here as this is the only sync point where we
295 know (without checking or additional sync-ing) that prune tasks in
296 in the current stream have completed (having just blocking-waited
297 for the force D2H). */
298 countPruneKernelTime(timers, timings, iLocality);
300 /* only count atdat and pair-list H2D at pair-search step */
301 if (timers->interaction[iLocality].didPairlistH2D)
303 /* atdat transfer timing (add only once, at local F wait) */
304 if (atomLocality == AtomLocality::Local)
307 timings->pl_h2d_t += timers->atdat.getLastRangeTime();
310 timings->pl_h2d_t += timers->interaction[iLocality].pl_h2d.getLastRangeTime();
312 /* Clear the timing flag for the next step */
313 timers->interaction[iLocality].didPairlistH2D = false;
317 //TODO: move into shared source file with gmx_compile_cpp_as_cuda
318 //NOLINTNEXTLINE(misc-definitions-in-headers)
319 bool gpu_try_finish_task(gmx_nbnxn_gpu_t *nb,
321 const AtomLocality aloc,
322 const bool haveOtherWork,
326 GpuTaskCompletion completionKind)
328 GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
330 /* determine interaction locality from atom locality */
331 const InteractionLocality iLocality = gpuAtomToInteractionLocality(aloc);
333 // We skip when during the non-local phase there was actually no work to do.
334 // This is consistent with nbnxn_gpu_launch_kernel.
335 if (haveOtherWork || !canSkipWork(*nb, iLocality))
337 // Query the state of the GPU stream and return early if we're not done
338 if (completionKind == GpuTaskCompletion::Check)
340 if (!haveStreamTasksCompleted(nb->stream[iLocality]))
342 // Early return to skip the steps below that we have to do only
343 // after the NB task completed
349 gpuStreamSynchronize(nb->stream[iLocality]);
352 bool calcEner = (flags & GMX_FORCE_ENERGY) != 0;
353 bool calcFshift = (flags & GMX_FORCE_VIRIAL) != 0;
355 gpu_accumulate_timings(nb->timings, nb->timers, nb->plist[iLocality], aloc, calcEner,
358 gpu_reduce_staged_outputs(nb->nbst, iLocality, calcEner, calcFshift, e_lj, e_el, fshift);
361 /* Always reset both pruning flags (doesn't hurt doing it even when timing is off). */
362 nb->timers->interaction[iLocality].didPrune = nb->timers->interaction[iLocality].didRollingPrune = false;
364 /* Turn off initial list pruning (doesn't hurt if this is not pair-search step). */
365 nb->plist[iLocality]->haveFreshList = false;
371 * Wait for the asynchronously launched nonbonded tasks and data
372 * transfers to finish.
374 * Also does timing accounting and reduction of the internal staging buffers.
375 * As this is called at the end of the step, it also resets the pair list and
378 * \param[in] nb The nonbonded data GPU structure
379 * \param[in] flags Force flags
380 * \param[in] aloc Atom locality identifier
381 * \param[in] haveOtherWork Tells whether there is other work than non-bonded work in the nbnxn stream(s)
382 * \param[out] e_lj Pointer to the LJ energy output to accumulate into
383 * \param[out] e_el Pointer to the electrostatics energy output to accumulate into
384 * \param[out] fshift Pointer to the shift force buffer to accumulate into
386 //NOLINTNEXTLINE(misc-definitions-in-headers) TODO: move into source file
387 void gpu_wait_finish_task(gmx_nbnxn_gpu_t *nb,
395 gpu_try_finish_task(nb, flags, aloc, haveOtherWork, e_lj, e_el, fshift,
396 GpuTaskCompletion::Wait);