Modernize wallcycle counting
[alexxy/gromacs.git] / src / gromacs / timing / wallcycle.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2008, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2017,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifndef GMX_TIMING_WALLCYCLE_H
39 #define GMX_TIMING_WALLCYCLE_H
40
41 /* NOTE: None of the routines here are safe to call within an OpenMP
42  * region */
43
44 #include <stdio.h>
45
46 #include <array>
47 #include <memory>
48 #include <vector>
49
50 #include "gromacs/timing/cyclecounter.h"
51 #include "gromacs/utility/basedefinitions.h"
52 #include "gromacs/utility/enumerationhelpers.h"
53
54 #ifndef DEBUG_WCYCLE
55 /*! \brief Enables consistency checking for the counters.
56  *
57  * If the macro is set to 1, code checks if you stop a counter different from the last
58  * one that was opened and if you do nest too deep.
59  */
60 #    define DEBUG_WCYCLE 0
61 #endif
62
63 struct t_commrec;
64
65 #ifndef DEBUG_WCYCLE
66 /*! \brief Enables consistency checking for the counters.
67  *
68  * If the macro is set to 1, code checks if you stop a counter different from the last
69  * one that was opened and if you do nest too deep.
70  */
71 #    define DEBUG_WCYCLE 0
72 #endif
73
74 enum class WallCycleCounter : int
75 {
76     Run,
77     Step,
78     PpDuringPme,
79     Domdec,
80     DDCommLoad,
81     DDCommBound,
82     VsiteConstr,
83     PpPmeSendX,
84     NS,
85     LaunchGpu,
86     MoveX,
87     Force,
88     MoveF,
89     PmeMesh,
90     PmeRedistXF,
91     PmeSpread,
92     PmeGather,
93     PmeFft,
94     PmeFftComm,
95     LJPme,
96     PmeSolve,
97     PmeWaitComm,
98     PpPmeWaitRecvF,
99     WaitGpuPmeSpread,
100     PmeFftMixedMode,
101     PmeSolveMixedMode,
102     WaitGpuPmeGather,
103     WaitGpuBonded,
104     PmeGpuFReduction,
105     WaitGpuNbNL,
106     WaitGpuNbL,
107     WaitGpuStatePropagatorData,
108     NbXFBufOps,
109     VsiteSpread,
110     PullPot,
111     Awh,
112     Traj,
113     Update,
114     Constr,
115     MoveE,
116     Rot,
117     RotAdd,
118     Swap,
119     Imd,
120     Test,
121     Count
122 };
123
124 enum class WallCycleSubCounter : int
125 {
126     DDRedist,
127     DDGrid,
128     DDSetupComm,
129     DDMakeTop,
130     DDMakeConstr,
131     DDTopOther,
132     DDGpu,
133     NBSGridLocal,
134     NBSGridNonLocal,
135     NBSSearchLocal,
136     NBSSearchNonLocal,
137     Listed,
138     ListedFep,
139     Restraints,
140     ListedBufOps,
141     NonbondedPruning,
142     NonbondedKernel,
143     NonbondedClear,
144     NonbondedFep,
145     LaunchGpuNonBonded,
146     LaunchGpuBonded,
147     LaunchGpuPme,
148     LaunchStatePropagatorData,
149     EwaldCorrection,
150     NBXBufOps,
151     NBFBufOps,
152     ClearForceBuffer,
153     LaunchGpuNBXBufOps,
154     LaunchGpuNBFBufOps,
155     LaunchGpuMoveX,
156     LaunchGpuMoveF,
157     LaunchGpuUpdateConstrain,
158     Test,
159     Count
160 };
161
162 static constexpr int sc_numWallCycleCounters        = static_cast<int>(WallCycleCounter::Count);
163 static constexpr int sc_numWallCycleSubCounters     = static_cast<int>(WallCycleSubCounter::Count);
164 static constexpr int sc_numWallCycleCountersSquared = sc_numWallCycleCounters * sc_numWallCycleCounters;
165 static constexpr bool sc_useCycleSubcounters        = GMX_CYCLE_SUBCOUNTERS;
166
167 struct wallcc_t
168 {
169     int          n;
170     gmx_cycles_t c;
171     gmx_cycles_t start;
172 };
173
174 #if DEBUG_WCYCLE
175 static constexpr int c_MaxWallCycleDepth = 6;
176 #endif
177
178
179 struct gmx_wallcycle
180 {
181     gmx::EnumerationArray<WallCycleCounter, wallcc_t> wcc;
182     /* did we detect one or more invalid cycle counts */
183     bool haveInvalidCount;
184     /* variables for testing/debugging */
185     bool                  wc_barrier;
186     std::vector<wallcc_t> wcc_all;
187     int                   wc_depth;
188 #if DEBUG_WCYCLE
189     std::array<WallCycleCounter, c_MaxWallCycleDepth> counterlist;
190     int                                               count_depth;
191     bool                                              isMasterRank;
192 #endif
193     WallCycleCounter                                     ewc_prev;
194     gmx_cycles_t                                         cycle_prev;
195     int64_t                                              reset_counters;
196     const t_commrec*                                     cr;
197     gmx::EnumerationArray<WallCycleSubCounter, wallcc_t> wcsc;
198 };
199
200 //! Returns if cycle counting is supported
201 bool wallcycle_have_counter();
202
203 //! Returns the wall cycle structure.
204 std::unique_ptr<gmx_wallcycle> wallcycle_init(FILE* fplog, int resetstep, const t_commrec* cr);
205
206 //! Adds custom barrier for wallcycle counting.
207 void wallcycleBarrier(gmx_wallcycle* wc);
208
209 void wallcycle_sub_get(gmx_wallcycle* wc, WallCycleSubCounter ewcs, int* n, double* c);
210 /* Returns the cumulative count and sub cycle count for ewcs */
211
212 inline void wallcycle_all_start(gmx_wallcycle* wc, WallCycleCounter ewc, gmx_cycles_t cycle)
213 {
214     wc->ewc_prev   = ewc;
215     wc->cycle_prev = cycle;
216 }
217
218 inline void wallcycle_all_stop(gmx_wallcycle* wc, WallCycleCounter ewc, gmx_cycles_t cycle)
219 {
220     const int prev    = static_cast<int>(wc->ewc_prev);
221     const int current = static_cast<int>(ewc);
222     wc->wcc_all[prev * sc_numWallCycleCounters + current].n += 1;
223     wc->wcc_all[prev * sc_numWallCycleCounters + current].c += cycle - wc->cycle_prev;
224 }
225
226 //! Starts the cycle counter (and increases the call count)
227 inline void wallcycle_start(gmx_wallcycle* wc, WallCycleCounter ewc)
228 {
229     if (wc == nullptr)
230     {
231         return;
232     }
233
234     wallcycleBarrier(wc);
235
236 #if DEBUG_WCYCLE
237     debug_start_check(wc, ewc);
238 #endif
239     gmx_cycles_t cycle = gmx_cycles_read();
240     wc->wcc[ewc].start = cycle;
241     if (!wc->wcc_all.empty())
242     {
243         wc->wc_depth++;
244         if (ewc == WallCycleCounter::Run)
245         {
246             wallcycle_all_start(wc, ewc, cycle);
247         }
248         else if (wc->wc_depth == 3)
249         {
250             wallcycle_all_stop(wc, ewc, cycle);
251         }
252     }
253 }
254
255 //! Starts the cycle counter without increasing the call count
256 inline void wallcycle_start_nocount(gmx_wallcycle* wc, WallCycleCounter ewc)
257 {
258     if (wc == nullptr)
259     {
260         return;
261     }
262     wc->wcc[ewc].n++;
263 }
264
265 //! Stop the cycle count for ewc , returns the last cycle count
266 inline double wallcycle_stop(gmx_wallcycle* wc, WallCycleCounter ewc)
267 {
268     gmx_cycles_t cycle, last;
269
270     if (wc == nullptr)
271     {
272         return 0;
273     }
274
275     wallcycleBarrier(wc);
276
277 #if DEBUG_WCYCLE
278     debug_stop_check(wc, ewc);
279 #endif
280
281     /* When processes or threads migrate between cores, the cycle counting
282      * can get messed up if the cycle counter on different cores are not
283      * synchronized. When this happens we expect both large negative and
284      * positive cycle differences. We can detect negative cycle differences.
285      * Detecting too large positive counts if difficult, since count can be
286      * large, especially for ewcRUN. If we detect a negative count,
287      * we will not print the cycle accounting table.
288      */
289     cycle = gmx_cycles_read();
290     if (cycle >= wc->wcc[ewc].start)
291     {
292         last = cycle - wc->wcc[ewc].start;
293     }
294     else
295     {
296         last                 = 0;
297         wc->haveInvalidCount = true;
298     }
299     wc->wcc[ewc].c += last;
300     wc->wcc[ewc].n++;
301     if (!wc->wcc_all.empty())
302     {
303         wc->wc_depth--;
304         if (ewc == WallCycleCounter::Run)
305         {
306             wallcycle_all_stop(wc, ewc, cycle);
307         }
308         else if (wc->wc_depth == 2)
309         {
310             wallcycle_all_start(wc, ewc, cycle);
311         }
312     }
313
314     return last;
315 }
316
317 //! Only increment call count for ewc by one
318 inline void wallcycle_increment_event_count(gmx_wallcycle* wc, WallCycleCounter ewc)
319 {
320     if (wc == nullptr)
321     {
322         return;
323     }
324     wc->wcc[ewc].n++;
325 }
326
327 //! Returns the cumulative count and cycle count for ewc
328 void wallcycle_get(gmx_wallcycle* wc, WallCycleCounter ewc, int* n, double* c);
329
330 //! Resets all cycle counters to zero
331 void wallcycle_reset_all(gmx_wallcycle* wc);
332
333 //! Scale the cycle counts to reflect how many threads run for that number of cycles
334 void wallcycle_scale_by_num_threads(gmx_wallcycle* wc, bool isPmeRank, int nthreads_pp, int nthreads_pme);
335
336 //! Return reset_counters from wc struct
337 int64_t wcycle_get_reset_counters(gmx_wallcycle* wc);
338
339 //! Set reset_counters
340 void wcycle_set_reset_counters(gmx_wallcycle* wc, int64_t reset_counters);
341
342 //! Set the start sub cycle count for ewcs
343 inline void wallcycle_sub_start(gmx_wallcycle* wc, WallCycleSubCounter ewcs)
344 {
345     if (sc_useCycleSubcounters && wc != nullptr)
346     {
347         wc->wcsc[ewcs].start = gmx_cycles_read();
348     }
349 }
350
351 //! Set the start sub cycle count for ewcs without increasing the call count
352 inline void wallcycle_sub_start_nocount(gmx_wallcycle* wc, WallCycleSubCounter ewcs)
353 {
354     if (sc_useCycleSubcounters && wc != nullptr)
355     {
356         wc->wcsc[ewcs].start = gmx_cycles_read();
357     }
358 }
359
360 //! Stop the sub cycle count for ewcs
361 inline void wallcycle_sub_stop(gmx_wallcycle* wc, WallCycleSubCounter ewcs)
362 {
363     if (sc_useCycleSubcounters && wc != nullptr)
364     {
365         wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
366         wc->wcsc[ewcs].n++;
367     }
368 }
369
370 #endif