Implement alternating GPU wait
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_gpu_common.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017, 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_mdlib
41  */
42
43 #ifndef GMX_MDLIB_NBNXN_GPU_COMMON_H
44 #define GMX_MDLIB_NBNXN_GPU_COMMON_H
45
46 #include "config.h"
47
48 #include <string>
49
50 #if GMX_GPU == GMX_GPU_CUDA
51 #include "nbnxn_cuda/nbnxn_cuda_types.h"
52 #endif
53
54 #if GMX_GPU == GMX_GPU_OPENCL
55 #include "nbnxn_ocl/nbnxn_ocl_types.h"
56 #endif
57
58 #include "gromacs/gpu_utils/gpu_utils.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdlib/nbnxn_gpu_types.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/timing/gpu_timing.h"
63 #include "gromacs/utility/stringutil.h"
64
65 #include "nbnxn_gpu_common_utils.h"
66
67 /*! \brief Check that atom locality values are valid for the GPU module.
68  *
69  *  In the GPU module atom locality "all" is not supported, the local and
70  *  non-local ranges are treated separately.
71  *
72  *  \param[in] atomLocality atom locality specifier
73  */
74 static inline void validateGpuAtomLocality(int atomLocality)
75 {
76     std::string str = gmx::formatString("Invalid atom locality passed (%d); valid here is only "
77                                         "local (%d) or nonlocal (%d)", atomLocality, eatLocal, eatNonlocal);
78
79     GMX_ASSERT(LOCAL_OR_NONLOCAL_A(atomLocality), str.c_str());
80 }
81
82 /*! \brief Convert atom locality to interaction locality.
83  *
84  *  In the current implementation the this is straightforward conversion:
85  *  local to local, non-local to non-local.
86  *
87  *  \param[in] atomLocality Atom locality specifier
88  *  \returns                Interaction locality corresponding to the atom locality passed.
89  */
90 static inline int gpuAtomToInteractionLocality(int atomLocality)
91 {
92     validateGpuAtomLocality(atomLocality);
93
94     /* determine interaction locality from atom locality */
95     if (LOCAL_A(atomLocality))
96     {
97         return eintLocal;
98     }
99     else if (NONLOCAL_A(atomLocality))
100     {
101         return eintNonlocal;
102     }
103     else
104     {
105         // can't be reached
106         assert(false);
107         return -1;
108     }
109 }
110
111 /*! \brief Calculate atom range and return start index and length.
112  *
113  * \param[in] atomData Atom descriptor data structure
114  * \param[in] atomLocality Atom locality specifier
115  * \param[out] atomRangeBegin Starting index of the atom range in the atom data array.
116  * \param[out] atomRangeLen Atom range length in the atom data array.
117  */
118 template <typename AtomDataT>
119 static inline void getGpuAtomRange(const AtomDataT *atomData,
120                                    int              atomLocality,
121                                    int             &atomRangeBegin,
122                                    int             &atomRangeLen)
123 {
124     assert(atomData);
125     validateGpuAtomLocality(atomLocality);
126
127     /* calculate the atom data index range based on locality */
128     if (LOCAL_A(atomLocality))
129     {
130         atomRangeBegin  = 0;
131         atomRangeLen    = atomData->natoms_local;
132     }
133     else
134     {
135         atomRangeBegin  = atomData->natoms_local;
136         atomRangeLen    = atomData->natoms - atomData->natoms_local;
137     }
138 }
139
140
141 /*! \brief Count pruning kernel time if either kernel has been triggered
142  *
143  *  We do the accounting for either of the two pruning kernel flavors:
144  *   - 1st pass prune: ran during the current step (prior to the force kernel);
145  *   - rolling prune:  ran at the end of the previous step (prior to the current step H2D xq);
146  *
147  * Note that the resetting of cu_timers_t::didPrune and cu_timers_t::didRollingPrune should happen
148  * after calling this function.
149  *
150  * \param[in] timers   structs with GPU timer objects
151  * \param[inout] timings  GPU task timing data
152  * \param[in] iloc        interaction locality
153  */
154 template <typename GpuTimers>
155 static void countPruneKernelTime(GpuTimers                 *timers,
156                                  gmx_wallclock_gpu_nbnxn_t *timings,
157                                  const int                  iloc)
158 {
159     // We might have not done any pruning (e.g. if we skipped with empty domains).
160     if (!timers->didPrune[iloc] && !timers->didRollingPrune[iloc])
161     {
162         return;
163     }
164
165     if (timers->didPrune[iloc])
166     {
167         timings->pruneTime.c++;
168         timings->pruneTime.t += timers->prune_k[iloc].getLastRangeTime();
169     }
170
171     if (timers->didRollingPrune[iloc])
172     {
173         timings->dynamicPruneTime.c++;
174         timings->dynamicPruneTime.t += timers->rollingPrune_k[iloc].getLastRangeTime();
175     }
176 }
177
178 /*! \brief Reduce data staged internally in the nbnxn module.
179  *
180  * Shift forces and electrostatic/LJ energies copied from the GPU into
181  * a module-internal staging area are immediately reduced (CPU-side buffers passed)
182  * after having waited for the transfers' completion.
183  *
184  * Note that this function should always be called after the transfers into the
185  * staging buffers has completed.
186  *
187  * \tparam     StagingData    Type of staging data
188  * \param[in]  nbst           Nonbonded staging data
189  * \param[in]  iLocality      Interaction locality specifier
190  * \param[in]  reduceEnergies True if energy reduction should be done
191  * \param[in]  reduceFshift   True if shift force reduction should be done
192  * \param[out] e_lj           Variable to accumulate LJ energy into
193  * \param[out] e_el           Variable to accumulate electrostatic energy into
194  * \param[out] fshift         Pointer to the array of shift forces to accumulate into
195  */
196 template <typename StagingData>
197 static inline void nbnxn_gpu_reduce_staged_outputs(const StagingData &nbst,
198                                                    int                iLocality,
199                                                    bool               reduceEnergies,
200                                                    bool               reduceFshift,
201                                                    real              *e_lj,
202                                                    real              *e_el,
203                                                    rvec              *fshift)
204 {
205     /* add up energies and shift forces (only once at local F wait) */
206     if (LOCAL_I(iLocality))
207     {
208         if (reduceEnergies)
209         {
210             *e_lj += *nbst.e_lj;
211             *e_el += *nbst.e_el;
212         }
213
214         if (reduceFshift)
215         {
216             for (int i = 0; i < SHIFTS; i++)
217             {
218                 rvec_inc(fshift[i], nbst.fshift[i]);
219             }
220         }
221     }
222 }
223
224 /*! \brief Do the per-step timing accounting of the nonbonded tasks.
225  *
226  *  Does timing accumulation and call-count increments for the nonbonded kernels.
227  *  Note that this function should be called after the current step's nonbonded
228  *  nonbonded tasks have completed with the exception of the rolling pruning kernels
229  *  that are accounted for during the following step.
230  *
231  * NOTE: if timing with multiple GPUs (streams) becomes possible, the
232  *      counters could end up being inconsistent due to not being incremented
233  *      on some of the node when this is skipped on empty local domains!
234  *
235  * \tparam     GpuTimers         GPU timers type
236  * \tparam     GpuPairlist       Pair list type
237  * \param[out] timings           Pointer to the NB GPU timings data
238  * \param[in]  timers            Pointer to GPU timers data
239  * \param[in]  plist             Pointer to the pair list data
240  * \param[in]  atomLocality      Atom locality specifier
241  * \param[in]  didEnergyKernels  True if energy kernels have been called in the current step
242  * \param[in]  doTiming          True if timing is enabled.
243  *
244  */
245 template <typename GpuTimers, typename GpuPairlist>
246 static inline void nbnxn_gpu_accumulate_timings(gmx_wallclock_gpu_nbnxn_t *timings,
247                                                 GpuTimers                 *timers,
248                                                 const GpuPairlist         *plist,
249                                                 int                        atomLocality,
250                                                 bool                       didEnergyKernels,
251                                                 bool                       doTiming)
252 {
253     /* timing data accumulation */
254     if (!doTiming)
255     {
256         return;
257     }
258
259     /* determine interaction locality from atom locality */
260     int iLocality = gpuAtomToInteractionLocality(atomLocality);
261
262     /* only increase counter once (at local F wait) */
263     if (LOCAL_I(iLocality))
264     {
265         timings->nb_c++;
266         timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].c += 1;
267     }
268
269     /* kernel timings */
270     timings->ktime[plist->haveFreshList ? 1 : 0][didEnergyKernels ? 1 : 0].t +=
271         timers->nb_k[iLocality].getLastRangeTime();
272
273     /* X/q H2D and F D2H timings */
274     timings->nb_h2d_t += timers->nb_h2d[iLocality].getLastRangeTime();
275     timings->nb_d2h_t += timers->nb_d2h[iLocality].getLastRangeTime();
276
277     /* Count the pruning kernel times for both cases:1st pass (at search step)
278        and rolling pruning (if called at the previous step).
279        We do the accounting here as this is the only sync point where we
280        know (without checking or additional sync-ing) that prune tasks in
281        in the current stream have completed (having just blocking-waited
282        for the force D2H). */
283     countPruneKernelTime(timers, timings, iLocality);
284
285     /* only count atdat and pair-list H2D at pair-search step */
286     if (timers->didPairlistH2D[iLocality])
287     {
288         /* atdat transfer timing (add only once, at local F wait) */
289         if (LOCAL_A(atomLocality))
290         {
291             timings->pl_h2d_c++;
292             timings->pl_h2d_t += timers->atdat.getLastRangeTime();
293         }
294
295         timings->pl_h2d_t += timers->pl_h2d[iLocality].getLastRangeTime();
296
297         /* Clear the timing flag for the next step */
298         timers->didPairlistH2D[iLocality] = false;
299     }
300 }
301
302 bool nbnxn_gpu_try_finish_task(gmx_nbnxn_gpu_t  *nb,
303                                int               flags,
304                                int               aloc,
305                                real             *e_lj,
306                                real             *e_el,
307                                rvec             *fshift,
308                                GpuTaskCompletion completionKind)
309 {
310     /* determine interaction locality from atom locality */
311     int iLocality = gpuAtomToInteractionLocality(aloc);
312
313     //  We skip when during the non-local phase there was actually no work to do.
314     //  This is consistent with nbnxn_gpu_launch_kernel.
315     if (!canSkipWork(nb, iLocality))
316     {
317         // Query the state of the GPU stream and return early if we're not done
318         if (completionKind == GpuTaskCompletion::Check)
319         {
320             if (!haveStreamTasksCompleted(nb->stream[iLocality]))
321             {
322                 // Early return to skip the steps below that we have to do only
323                 // after the NB task completed
324                 return false;
325             }
326         }
327         else
328         {
329             gpuStreamSynchronize(nb->stream[iLocality]);
330         }
331
332         bool calcEner   = flags & GMX_FORCE_ENERGY;
333         bool calcFshift = flags & GMX_FORCE_VIRIAL;
334
335         nbnxn_gpu_accumulate_timings(nb->timings, nb->timers, nb->plist[iLocality], aloc, calcEner, nb->bDoTime);
336
337         nbnxn_gpu_reduce_staged_outputs(nb->nbst, iLocality, calcEner, calcFshift, e_lj, e_el, fshift);
338     }
339
340     /* Always reset both pruning flags (doesn't hurt doing it even when timing is off). */
341     nb->timers->didPrune[iLocality] = nb->timers->didRollingPrune[iLocality] = false;
342
343     /* Turn off initial list pruning (doesn't hurt if this is not pair-search step). */
344     nb->plist[iLocality]->haveFreshList = false;
345
346     return true;
347 }
348
349 /*! \brief
350  * Wait for the asynchronously launched nonbonded tasks and data
351  * transfers to finish.
352  *
353  * Also does timing accounting and reduction of the internal staging buffers.
354  * As this is called at the end of the step, it also resets the pair list and
355  * pruning flags.
356  *
357  * \param[in] nb The nonbonded data GPU structure
358  * \param[in] flags Force flags
359  * \param[in] aloc Atom locality identifier
360  * \param[out] e_lj Pointer to the LJ energy output to accumulate into
361  * \param[out] e_el Pointer to the electrostatics energy output to accumulate into
362  * \param[out] fshift Pointer to the shift force buffer to accumulate into
363  */
364 void nbnxn_gpu_wait_finish_task(gmx_nbnxn_gpu_t *nb,
365                                 int              flags,
366                                 int              aloc,
367                                 real            *e_lj,
368                                 real            *e_el,
369                                 rvec            *fshift)
370 {
371     nbnxn_gpu_try_finish_task(nb, flags, aloc, e_lj, e_el, fshift,
372                               GpuTaskCompletion::Wait);
373 }
374
375 #endif