Replace rounding using int(x+0.5) with roundToInt
[alexxy/gromacs.git] / src / gromacs / timing / wallcycle.cpp
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,2016,2017,2018, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "wallcycle.h"
40
41 #include "config.h"
42
43 #include <cstdlib>
44
45 #include <array>
46 #include <vector>
47
48 #include "gromacs/math/functions.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/timing/cyclecounter.h"
51 #include "gromacs/timing/gpu_timing.h"
52 #include "gromacs/timing/wallcyclereporting.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/gmxmpi.h"
56 #include "gromacs/utility/logger.h"
57 #include "gromacs/utility/smalloc.h"
58 #include "gromacs/utility/snprintf.h"
59
60 static const bool useCycleSubcounters = GMX_CYCLE_SUBCOUNTERS;
61
62 /* DEBUG_WCYCLE adds consistency checking for the counters.
63  * It checks if you stop a counter different from the last
64  * one that was opened and if you do nest too deep.
65  */
66 /* #define DEBUG_WCYCLE */
67
68 #ifdef DEBUG_WCYCLE
69 #include "gromacs/utility/fatalerror.h"
70 #endif
71
72 typedef struct
73 {
74     int          n;
75     gmx_cycles_t c;
76     gmx_cycles_t start;
77 } wallcc_t;
78
79 struct gmx_wallcycle
80 {
81     wallcc_t        *wcc;
82     /* did we detect one or more invalid cycle counts */
83     gmx_bool         haveInvalidCount;
84     /* variables for testing/debugging */
85     gmx_bool         wc_barrier;
86     wallcc_t        *wcc_all;
87     int              wc_depth;
88 #ifdef DEBUG_WCYCLE
89 #define DEPTH_MAX 6
90     int               counterlist[DEPTH_MAX];
91     int               count_depth;
92 #endif
93     int               ewc_prev;
94     gmx_cycles_t      cycle_prev;
95     int64_t           reset_counters;
96 #if GMX_MPI
97     MPI_Comm          mpi_comm_mygroup;
98 #endif
99     wallcc_t         *wcsc;
100 };
101
102 /* Each name should not exceed 19 printing characters
103    (ie. terminating null can be twentieth) */
104 static const char *wcn[ewcNR] =
105 {
106     "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
107     "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
108     "Comm. coord.", "Force", "Wait + Comm. F", "PME mesh",
109     "PME redist. X/F", "PME spread", "PME gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve LJ", "PME solve Elec",
110     "PME wait for PP", "Wait + Recv. PME F",
111     "Wait PME GPU spread", "PME 3D-FFT", "PME solve", /* the strings for FFT/solve are repeated here for mixed mode counters */
112     "Wait PME GPU gather", "Reduce GPU PME F",
113     "Wait GPU NB nonloc.", "Wait GPU NB local", "NB X/F buffer ops.",
114     "Vsite spread", "COM pull force", "AWH",
115     "Write traj.", "Update", "Constraints", "Comm. energies",
116     "Enforced rotation", "Add rot. forces", "Position swapping", "IMD", "Test"
117 };
118
119 static const char *wcsn[ewcsNR] =
120 {
121     "DD redist.", "DD NS grid + sort", "DD setup comm.",
122     "DD make top.", "DD make constr.", "DD top. other",
123     "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
124     "Bonded F",
125     "Bonded-FEP F",
126     "Restraints F",
127     "Listed buffer ops.",
128     "Nonbonded pruning",
129     "Nonbonded F",
130     "Launch NB GPU tasks",
131     "Launch PME GPU tasks",
132     "Ewald F correction",
133     "NB X buffer ops.",
134     "NB F buffer ops.",
135 };
136
137 /* PME GPU timing events' names - correspond to the enum in the gpu_timing.h */
138 static const char *PMEStageNames[] =
139 {
140     "PME spline",
141     "PME spread",
142     "PME spline + spread",
143     "PME 3D-FFT r2c",
144     "PME solve",
145     "PME 3D-FFT c2r",
146     "PME gather",
147 };
148
149 gmx_bool wallcycle_have_counter()
150 {
151     return gmx_cycles_have_counter();
152 }
153
154 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec gmx_unused *cr)
155 {
156     gmx_wallcycle_t wc;
157
158
159     if (!wallcycle_have_counter())
160     {
161         return nullptr;
162     }
163
164     snew(wc, 1);
165
166     wc->haveInvalidCount    = FALSE;
167     wc->wc_barrier          = FALSE;
168     wc->wcc_all             = nullptr;
169     wc->wc_depth            = 0;
170     wc->ewc_prev            = -1;
171     wc->reset_counters      = resetstep;
172
173 #if GMX_MPI
174     if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != nullptr)
175     {
176         if (fplog)
177         {
178             fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
179         }
180         wc->wc_barrier       = TRUE;
181         wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
182     }
183 #endif
184
185     snew(wc->wcc, ewcNR);
186     if (getenv("GMX_CYCLE_ALL") != nullptr)
187     {
188         if (fplog)
189         {
190             fprintf(fplog, "\nWill time all the code during the run\n\n");
191         }
192         snew(wc->wcc_all, ewcNR*ewcNR);
193     }
194
195     if (useCycleSubcounters)
196     {
197         snew(wc->wcsc, ewcsNR);
198     }
199
200 #ifdef DEBUG_WCYCLE
201     wc->count_depth = 0;
202 #endif
203
204     return wc;
205 }
206
207 /* TODO: Should be called from finish_run() or runner()
208    void wallcycle_destroy(gmx_wallcycle_t wc)
209    {
210     if (wc == nullptr)
211     {
212         return;
213     }
214
215     if (wc->wcc != nullptr)
216     {
217         sfree(wc->wcc);
218     }
219     if (wc->wcc_all != nullptr)
220     {
221         sfree(wc->wcc_all);
222     }
223     if (wc->wcsc != nullptr)
224     {
225         sfree(wc->wcsc);
226     }
227     sfree(wc);
228    }
229  */
230
231 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
232 {
233     wc->ewc_prev   = ewc;
234     wc->cycle_prev = cycle;
235 }
236
237 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
238 {
239     wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
240     wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
241 }
242
243
244 #ifdef DEBUG_WCYCLE
245 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
246 {
247     /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
248
249     if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
250     {
251         gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
252                   wc->count_depth);
253     }
254     wc->counterlist[wc->count_depth] = ewc;
255     wc->count_depth++;
256 }
257
258 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
259 {
260     wc->count_depth--;
261
262     /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
263
264     if (wc->count_depth < 0)
265     {
266         gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
267     }
268     if (wc->counterlist[wc->count_depth] != ewc)
269     {
270         gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
271                   wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
272     }
273 }
274 #endif
275
276 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
277 {
278     gmx_cycles_t cycle;
279
280     if (wc == nullptr)
281     {
282         return;
283     }
284
285 #if GMX_MPI
286     if (wc->wc_barrier)
287     {
288         MPI_Barrier(wc->mpi_comm_mygroup);
289     }
290 #endif
291
292 #ifdef DEBUG_WCYCLE
293     debug_start_check(wc, ewc);
294 #endif
295
296     cycle              = gmx_cycles_read();
297     wc->wcc[ewc].start = cycle;
298     if (wc->wcc_all != nullptr)
299     {
300         wc->wc_depth++;
301         if (ewc == ewcRUN)
302         {
303             wallcycle_all_start(wc, ewc, cycle);
304         }
305         else if (wc->wc_depth == 3)
306         {
307             wallcycle_all_stop(wc, ewc, cycle);
308         }
309     }
310 }
311
312 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
313 {
314     if (wc == nullptr)
315     {
316         return;
317     }
318
319     wallcycle_start(wc, ewc);
320     wc->wcc[ewc].n--;
321 }
322
323 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
324 {
325     gmx_cycles_t cycle, last;
326
327     if (wc == nullptr)
328     {
329         return 0;
330     }
331
332 #if GMX_MPI
333     if (wc->wc_barrier)
334     {
335         MPI_Barrier(wc->mpi_comm_mygroup);
336     }
337 #endif
338
339 #ifdef DEBUG_WCYCLE
340     debug_stop_check(wc, ewc);
341 #endif
342
343     /* When processes or threads migrate between cores, the cycle counting
344      * can get messed up if the cycle counter on different cores are not
345      * synchronized. When this happens we expect both large negative and
346      * positive cycle differences. We can detect negative cycle differences.
347      * Detecting too large positive counts if difficult, since count can be
348      * large, especially for ewcRUN. If we detect a negative count,
349      * we will not print the cycle accounting table.
350      */
351     cycle                    = gmx_cycles_read();
352     if (cycle >= wc->wcc[ewc].start)
353     {
354         last                 = cycle - wc->wcc[ewc].start;
355     }
356     else
357     {
358         last                 = 0;
359         wc->haveInvalidCount = TRUE;
360     }
361     wc->wcc[ewc].c          += last;
362     wc->wcc[ewc].n++;
363     if (wc->wcc_all)
364     {
365         wc->wc_depth--;
366         if (ewc == ewcRUN)
367         {
368             wallcycle_all_stop(wc, ewc, cycle);
369         }
370         else if (wc->wc_depth == 2)
371         {
372             wallcycle_all_start(wc, ewc, cycle);
373         }
374     }
375
376     return last;
377 }
378
379 void wallcycle_get(gmx_wallcycle_t wc, int ewc, int *n, double *c)
380 {
381     *n = wc->wcc[ewc].n;
382     *c = static_cast<double>(wc->wcc[ewc].c);
383 }
384
385 void wallcycle_reset_all(gmx_wallcycle_t wc)
386 {
387     int i;
388
389     if (wc == nullptr)
390     {
391         return;
392     }
393
394     for (i = 0; i < ewcNR; i++)
395     {
396         wc->wcc[i].n = 0;
397         wc->wcc[i].c = 0;
398     }
399     wc->haveInvalidCount = FALSE;
400
401     if (wc->wcc_all)
402     {
403         for (i = 0; i < ewcNR*ewcNR; i++)
404         {
405             wc->wcc_all[i].n = 0;
406             wc->wcc_all[i].c = 0;
407         }
408     }
409     if (wc->wcsc)
410     {
411         for (i = 0; i < ewcsNR; i++)
412         {
413             wc->wcsc[i].n = 0;
414             wc->wcsc[i].c = 0;
415         }
416     }
417 }
418
419 static gmx_bool is_pme_counter(int ewc)
420 {
421     return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
422 }
423
424 static gmx_bool is_pme_subcounter(int ewc)
425 {
426     return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
427 }
428
429 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
430 static void subtract_cycles(wallcc_t *wcc, int ewc_main, int ewc_sub)
431 {
432     if (wcc[ewc_sub].n > 0)
433     {
434         if (wcc[ewc_main].c >= wcc[ewc_sub].c)
435         {
436             wcc[ewc_main].c -= wcc[ewc_sub].c;
437         }
438         else
439         {
440             /* Something is wrong with the cycle counting */
441             wcc[ewc_main].c  = 0;
442         }
443     }
444 }
445
446 void wallcycle_scale_by_num_threads(gmx_wallcycle_t wc, bool isPmeRank, int nthreads_pp, int nthreads_pme)
447 {
448     if (wc == nullptr)
449     {
450         return;
451     }
452
453     for (int i = 0; i < ewcNR; i++)
454     {
455         if (is_pme_counter(i) || (i == ewcRUN && isPmeRank))
456         {
457             wc->wcc[i].c *= nthreads_pme;
458
459             if (wc->wcc_all)
460             {
461                 for (int j = 0; j < ewcNR; j++)
462                 {
463                     wc->wcc_all[i*ewcNR+j].c *= nthreads_pme;
464                 }
465             }
466         }
467         else
468         {
469             wc->wcc[i].c *= nthreads_pp;
470
471             if (wc->wcc_all)
472             {
473                 for (int j = 0; j < ewcNR; j++)
474                 {
475                     wc->wcc_all[i*ewcNR+j].c *= nthreads_pp;
476                 }
477             }
478         }
479     }
480     if (useCycleSubcounters && wc->wcsc && !isPmeRank)
481     {
482         for (int i = 0; i < ewcsNR; i++)
483         {
484             wc->wcsc[i].c *= nthreads_pp;
485         }
486     }
487 }
488
489 /* TODO Make an object for this function to return, containing some
490  * vectors of something like wallcc_t for the summed wcc, wcc_all and
491  * wcsc, AND the original wcc for rank 0.
492  *
493  * The GPU timing is reported only for rank 0, so we want to preserve
494  * the original wcycle on that rank. Rank 0 also reports the global
495  * counts before that, so needs something to contain the global data
496  * without over-writing the rank-0 data. The current implementation
497  * uses cycles_sum to manage this, which works OK now because wcsc and
498  * wcc_all are unused by the GPU reporting, but it is not satisfactory
499  * for the future. Also, there's no need for MPI_Allreduce, since
500  * only MASTERRANK uses any of the results. */
501 WallcycleCounts wallcycle_sum(const t_commrec *cr, gmx_wallcycle_t wc)
502 {
503     WallcycleCounts cycles_sum;
504     wallcc_t       *wcc;
505     double          cycles[ewcNR+ewcsNR];
506 #if GMX_MPI
507     double          cycles_n[ewcNR+ewcsNR+1];
508 #endif
509     int             i;
510     int             nsum;
511
512     if (wc == nullptr)
513     {
514         /* Default construction of std::array of non-class T can leave
515            the values indeterminate, just like a C array, and icc
516            warns about it. */
517         cycles_sum.fill(0);
518         return cycles_sum;
519     }
520
521     wcc = wc->wcc;
522
523     subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
524     subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
525
526     subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
527
528     if (cr->npmenodes == 0)
529     {
530         /* All nodes do PME (or no PME at all) */
531         subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
532     }
533     else
534     {
535         /* The are PME-only nodes */
536         if (wcc[ewcPMEMESH].n > 0)
537         {
538             /* This must be a PME only node, calculate the Wait + Comm. time */
539             GMX_ASSERT(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c, "Total run ticks must be greater than PME-only ticks");
540             wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
541         }
542     }
543
544     /* Store the cycles in a double buffer for summing */
545     for (i = 0; i < ewcNR; i++)
546     {
547 #if GMX_MPI
548         cycles_n[i] = static_cast<double>(wcc[i].n);
549 #endif
550         cycles[i]   = static_cast<double>(wcc[i].c);
551     }
552     nsum = ewcNR;
553     if (wc->wcsc)
554     {
555         for (i = 0; i < ewcsNR; i++)
556         {
557 #if GMX_MPI
558             cycles_n[ewcNR+i] = static_cast<double>(wc->wcsc[i].n);
559 #endif
560             cycles[ewcNR+i]   = static_cast<double>(wc->wcsc[i].c);
561         }
562         nsum += ewcsNR;
563     }
564
565 #if GMX_MPI
566     if (cr->nnodes > 1)
567     {
568         double buf[ewcNR+ewcsNR+1];
569
570         // TODO this code is used only at the end of the run, so we
571         // can just do a simple reduce of haveInvalidCount in
572         // wallcycle_print, and avoid bugs
573         cycles_n[nsum] = (wc->haveInvalidCount ? 1 : 0);
574         // TODO Use MPI_Reduce
575         MPI_Allreduce(cycles_n, buf, nsum + 1, MPI_DOUBLE, MPI_MAX,
576                       cr->mpi_comm_mysim);
577         for (i = 0; i < ewcNR; i++)
578         {
579             wcc[i].n = gmx::roundToInt(buf[i]);
580         }
581         wc->haveInvalidCount = (buf[nsum] > 0);
582         if (wc->wcsc)
583         {
584             for (i = 0; i < ewcsNR; i++)
585             {
586                 wc->wcsc[i].n = gmx::roundToInt(buf[ewcNR+i]);
587             }
588         }
589
590         // TODO Use MPI_Reduce
591         MPI_Allreduce(cycles, cycles_sum.data(), nsum, MPI_DOUBLE, MPI_SUM,
592                       cr->mpi_comm_mysim);
593
594         if (wc->wcc_all != nullptr)
595         {
596             double *buf_all, *cyc_all;
597
598             snew(cyc_all, ewcNR*ewcNR);
599             snew(buf_all, ewcNR*ewcNR);
600             for (i = 0; i < ewcNR*ewcNR; i++)
601             {
602                 cyc_all[i] = wc->wcc_all[i].c;
603             }
604             // TODO Use MPI_Reduce
605             MPI_Allreduce(cyc_all, buf_all, ewcNR*ewcNR, MPI_DOUBLE, MPI_SUM,
606                           cr->mpi_comm_mysim);
607             for (i = 0; i < ewcNR*ewcNR; i++)
608             {
609                 wc->wcc_all[i].c = static_cast<gmx_cycles_t>(buf_all[i]);
610             }
611             sfree(buf_all);
612             sfree(cyc_all);
613         }
614     }
615     else
616 #endif
617     {
618         for (i = 0; i < nsum; i++)
619         {
620             cycles_sum[i] = cycles[i];
621         }
622     }
623
624     return cycles_sum;
625 }
626
627 static void print_cycles(FILE *fplog, double c2t, const char *name,
628                          int nnodes, int nthreads,
629                          int ncalls, double c_sum, double tot)
630 {
631     char   nnodes_str[STRLEN];
632     char   nthreads_str[STRLEN];
633     char   ncalls_str[STRLEN];
634     double wallt;
635     double percentage = (tot > 0.) ? (100. * c_sum / tot) : 0.;
636
637     if (c_sum > 0)
638     {
639         if (ncalls > 0)
640         {
641             snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
642             if (nnodes < 0)
643             {
644                 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
645             }
646             else
647             {
648                 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
649             }
650             if (nthreads < 0)
651             {
652                 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
653             }
654             else
655             {
656                 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
657             }
658         }
659         else
660         {
661             nnodes_str[0]   = 0;
662             nthreads_str[0] = 0;
663             ncalls_str[0]   = 0;
664         }
665         /* Convert the cycle count to wallclock time for this task */
666         wallt = c_sum*c2t;
667
668         fprintf(fplog, " %-19.19s %4s %4s %10s  %10.3f %14.3f %5.1f\n",
669                 name, nnodes_str, nthreads_str, ncalls_str, wallt,
670                 c_sum*1e-9, percentage);
671     }
672 }
673
674 static void print_gputimes(FILE *fplog, const char *name,
675                            int n, double t, double tot_t)
676 {
677     char num[11];
678     char avg_perf[11];
679
680     if (n > 0)
681     {
682         snprintf(num, sizeof(num), "%10d", n);
683         snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
684     }
685     else
686     {
687         sprintf(num, "          ");
688         sprintf(avg_perf, "          ");
689     }
690     if (t != tot_t && tot_t > 0)
691     {
692         fprintf(fplog, " %-29s %10s%12.3f   %s   %5.1f\n",
693                 name, num, t/1000, avg_perf, 100 * t/tot_t);
694     }
695     else
696     {
697         fprintf(fplog, " %-29s %10s%12.3f   %s   %5.1f\n",
698                 name, "", t/1000, avg_perf, 100.0);
699     }
700 }
701
702 static void print_header(FILE *fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
703 {
704     int nrank_tot = nrank_pp + nrank_pme;
705     if (0 == nrank_pme)
706     {
707         fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
708         if (nth_pp > 1)
709         {
710             fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
711         }
712         /* Don't report doing PP+PME, because we can't tell here if
713          * this is RF, etc. */
714     }
715     else
716     {
717         fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
718         if (nth_pp > 1)
719         {
720             fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
721         }
722         fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
723         if (nth_pme > 1)
724         {
725             fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
726         }
727     }
728
729     fprintf(fplog, "\n\n");
730     fprintf(fplog, " Computing:          Num   Num      Call    Wall time         Giga-Cycles\n");
731     fprintf(fplog, "                     Ranks Threads  Count      (s)         total sum    %%\n");
732 }
733
734
735 void wallcycle_print(FILE *fplog, const gmx::MDLogger &mdlog, int nnodes, int npme,
736                      int nth_pp, int nth_pme, double realtime,
737                      gmx_wallcycle_t wc, const WallcycleCounts &cyc_sum,
738                      const gmx_wallclock_gpu_nbnxn_t *gpu_nbnxn_t,
739                      const gmx_wallclock_gpu_pme_t *gpu_pme_t)
740 {
741     double      tot, tot_for_pp, tot_for_rest, tot_cpu_overlap, gpu_cpu_ratio;
742     double      c2t, c2t_pp, c2t_pme = 0;
743     int         i, j, npp, nth_tot;
744     char        buf[STRLEN];
745     const char *hline = "-----------------------------------------------------------------------------";
746
747     if (wc == nullptr)
748     {
749         return;
750     }
751
752     GMX_ASSERT(nth_pp > 0, "Number of particle-particle threads must be >0");
753     GMX_ASSERT(nth_pme > 0, "Number of PME threads must be >0");
754     GMX_ASSERT(nnodes > 0, "Number of nodes must be >0");
755     GMX_ASSERT(npme >= 0, "Number of PME nodes cannot be negative");
756     npp     = nnodes - npme;
757     /* npme is the number of PME-only ranks used, and we always do PP work */
758     GMX_ASSERT(npp > 0, "Number of particle-particle nodes must be >0");
759
760     nth_tot = npp*nth_pp + npme*nth_pme;
761
762     /* When using PME-only nodes, the next line is valid for both
763        PP-only and PME-only nodes because they started ewcRUN at the
764        same time. */
765     tot        = cyc_sum[ewcRUN];
766     tot_for_pp = 0;
767
768     if (tot <= 0.0)
769     {
770         /* TODO This is heavy handed, but until someone reworks the
771            code so that it is provably robust with respect to
772            non-positive values for all possible timer and cycle
773            counters, there is less value gained from printing whatever
774            timing data might still be sensible for some non-Jenkins
775            run, than is lost from diagnosing Jenkins FP exceptions on
776            runs about whose execution time we don't care. */
777         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
778                 "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a time accounting",
779                 tot);
780         return;
781     }
782
783     if (wc->haveInvalidCount)
784     {
785         GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: Detected invalid cycle counts, probably because threads moved between CPU cores that do not have synchronized cycle counters. Will not print the cycle accounting.");
786         return;
787     }
788
789
790     /* Conversion factor from cycles to seconds */
791     c2t     = realtime/tot;
792     c2t_pp  = c2t * nth_tot / static_cast<double>(npp*nth_pp);
793     if (npme > 0)
794     {
795         c2t_pme = c2t * nth_tot / static_cast<double>(npme*nth_pme);
796     }
797     else
798     {
799         c2t_pme = 0;
800     }
801
802     fprintf(fplog, "\n     R E A L   C Y C L E   A N D   T I M E   A C C O U N T I N G\n\n");
803
804     print_header(fplog, npp, nth_pp, npme, nth_pme);
805
806     fprintf(fplog, "%s\n", hline);
807     for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
808     {
809         if (is_pme_subcounter(i))
810         {
811             /* Do not count these at all */
812         }
813         else if (npme > 0 && is_pme_counter(i))
814         {
815             /* Print timing information for PME-only nodes, but add an
816              * asterisk so the reader of the table can know that the
817              * walltimes are not meant to add up. The asterisk still
818              * fits in the required maximum of 19 characters. */
819             char buffer[STRLEN];
820             snprintf(buffer, STRLEN, "%s *", wcn[i]);
821             print_cycles(fplog, c2t_pme, buffer,
822                          npme, nth_pme,
823                          wc->wcc[i].n, cyc_sum[i], tot);
824         }
825         else
826         {
827             /* Print timing information when it is for a PP or PP+PME
828                node */
829             print_cycles(fplog, c2t_pp, wcn[i],
830                          npp, nth_pp,
831                          wc->wcc[i].n, cyc_sum[i], tot);
832             tot_for_pp += cyc_sum[i];
833         }
834     }
835     if (wc->wcc_all != nullptr)
836     {
837         for (i = 0; i < ewcNR; i++)
838         {
839             for (j = 0; j < ewcNR; j++)
840             {
841                 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
842                 print_cycles(fplog, c2t_pp, buf,
843                              npp, nth_pp,
844                              wc->wcc_all[i*ewcNR+j].n,
845                              wc->wcc_all[i*ewcNR+j].c,
846                              tot);
847             }
848         }
849     }
850     tot_for_rest = tot * npp * nth_pp / static_cast<double>(nth_tot);
851     print_cycles(fplog, c2t_pp, "Rest",
852                  npp, nth_pp,
853                  -1, tot_for_rest - tot_for_pp, tot);
854     fprintf(fplog, "%s\n", hline);
855     print_cycles(fplog, c2t, "Total",
856                  npp, nth_pp,
857                  -1, tot, tot);
858     fprintf(fplog, "%s\n", hline);
859
860     if (npme > 0)
861     {
862         fprintf(fplog,
863                 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
864                 "    twice the total reported, but the cycle count total and %% are correct.\n"
865                 "%s\n", hline);
866     }
867
868     if (wc->wcc[ewcPMEMESH].n > 0)
869     {
870         // A workaround to not print breakdown when no subcounters were recorded.
871         // TODO: figure out and record PME GPU counters (what to do with the waiting ones?)
872         std::vector<int> validPmeSubcounterIndices;
873         for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
874         {
875             if (is_pme_subcounter(i) && wc->wcc[i].n > 0)
876             {
877                 validPmeSubcounterIndices.push_back(i);
878             }
879         }
880
881         if (!validPmeSubcounterIndices.empty())
882         {
883             fprintf(fplog, " Breakdown of PME mesh computation\n");
884             fprintf(fplog, "%s\n", hline);
885             for (auto i : validPmeSubcounterIndices)
886             {
887                 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i],
888                              npme > 0 ? npme : npp, nth_pme,
889                              wc->wcc[i].n, cyc_sum[i], tot);
890             }
891             fprintf(fplog, "%s\n", hline);
892         }
893     }
894
895     if (useCycleSubcounters && wc->wcsc)
896     {
897         fprintf(fplog, " Breakdown of PP computation\n");
898         fprintf(fplog, "%s\n", hline);
899         for (i = 0; i < ewcsNR; i++)
900         {
901             print_cycles(fplog, c2t_pp, wcsn[i],
902                          npp, nth_pp,
903                          wc->wcsc[i].n, cyc_sum[ewcNR+i], tot);
904         }
905         fprintf(fplog, "%s\n", hline);
906     }
907
908     /* print GPU timing summary */
909     double tot_gpu = 0.0;
910     if (gpu_pme_t)
911     {
912         for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
913         {
914             tot_gpu += gpu_pme_t->timing[k].t;
915         }
916     }
917     if (gpu_nbnxn_t)
918     {
919         const char *k_log_str[2][2] = {
920             {"Nonbonded F kernel", "Nonbonded F+ene k."},
921             {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
922         };
923         tot_gpu += gpu_nbnxn_t->pl_h2d_t + gpu_nbnxn_t->nb_h2d_t + gpu_nbnxn_t->nb_d2h_t;
924
925         /* add up the kernel timings */
926         for (i = 0; i < 2; i++)
927         {
928             for (j = 0; j < 2; j++)
929             {
930                 tot_gpu += gpu_nbnxn_t->ktime[i][j].t;
931             }
932         }
933         tot_gpu += gpu_nbnxn_t->pruneTime.t;
934
935         tot_cpu_overlap = wc->wcc[ewcFORCE].c;
936         if (wc->wcc[ewcPMEMESH].n > 0)
937         {
938             tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
939         }
940         tot_cpu_overlap *= realtime*1000/tot; /* convert s to ms */
941
942         fprintf(fplog, "\n GPU timings\n%s\n", hline);
943         fprintf(fplog, " Computing:                         Count  Wall t (s)      ms/step       %c\n", '%');
944         fprintf(fplog, "%s\n", hline);
945         print_gputimes(fplog, "Pair list H2D",
946                        gpu_nbnxn_t->pl_h2d_c, gpu_nbnxn_t->pl_h2d_t, tot_gpu);
947         print_gputimes(fplog, "X / q H2D",
948                        gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_h2d_t, tot_gpu);
949
950         for (i = 0; i < 2; i++)
951         {
952             for (j = 0; j < 2; j++)
953             {
954                 if (gpu_nbnxn_t->ktime[i][j].c)
955                 {
956                     print_gputimes(fplog, k_log_str[i][j],
957                                    gpu_nbnxn_t->ktime[i][j].c, gpu_nbnxn_t->ktime[i][j].t, tot_gpu);
958                 }
959             }
960         }
961         if (gpu_pme_t)
962         {
963             for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
964             {
965                 if (gpu_pme_t->timing[k].c)
966                 {
967                     print_gputimes(fplog, PMEStageNames[k],
968                                    gpu_pme_t->timing[k].c,
969                                    gpu_pme_t->timing[k].t,
970                                    tot_gpu);
971                 }
972             }
973         }
974         if (gpu_nbnxn_t->pruneTime.c)
975         {
976             print_gputimes(fplog, "Pruning kernel", gpu_nbnxn_t->pruneTime.c, gpu_nbnxn_t->pruneTime.t, tot_gpu);
977         }
978         print_gputimes(fplog, "F D2H",  gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_d2h_t, tot_gpu);
979         fprintf(fplog, "%s\n", hline);
980         print_gputimes(fplog, "Total ", gpu_nbnxn_t->nb_c, tot_gpu, tot_gpu);
981         fprintf(fplog, "%s\n", hline);
982         if (gpu_nbnxn_t->dynamicPruneTime.c)
983         {
984             /* We print the dynamic pruning kernel timings after a separator
985              * and avoid adding it to tot_gpu as this is not in the force
986              * overlap. We print the fraction as relative to the rest.
987              */
988             print_gputimes(fplog, "*Dynamic pruning", gpu_nbnxn_t->dynamicPruneTime.c, gpu_nbnxn_t->dynamicPruneTime.t, tot_gpu);
989             fprintf(fplog, "%s\n", hline);
990         }
991         gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
992         if (gpu_nbnxn_t->nb_c > 0 && wc->wcc[ewcFORCE].n > 0)
993         {
994             fprintf(fplog, "\nAverage per-step force GPU/CPU evaluation time ratio: %.3f ms/%.3f ms = %.3f\n",
995                     tot_gpu/gpu_nbnxn_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
996                     gpu_cpu_ratio);
997         }
998
999         /* only print notes related to CPU-GPU load balance with PME */
1000         if (wc->wcc[ewcPMEMESH].n > 0)
1001         {
1002             fprintf(fplog, "For optimal resource utilization this ratio should be close to 1\n");
1003
1004             /* print note if the imbalance is high with PME case in which
1005              * CPU-GPU load balancing is possible */
1006             if (gpu_cpu_ratio < 0.8 || gpu_cpu_ratio > 1.25)
1007             {
1008                 /* Only the sim master calls this function, so always print to stderr */
1009                 if (gpu_cpu_ratio < 0.8)
1010                 {
1011                     if (npp > 1)
1012                     {
1013                         /* The user could have used -notunepme,
1014                          * but we currently can't check that here.
1015                          */
1016                         GMX_LOG(mdlog.warning).asParagraph().appendText(
1017                                 "NOTE: The CPU has >25% more load than the GPU. This imbalance wastes\n"
1018                                 "      GPU resources. Maybe the domain decomposition limits the PME tuning.\n"
1019                                 "      In that case, try setting the DD grid manually (-dd) or lowering -dds.");
1020                     }
1021                     else
1022                     {
1023                         /* We should not end up here, unless the box is
1024                          * too small for increasing the cut-off for PME tuning.
1025                          */
1026                         GMX_LOG(mdlog.warning).asParagraph().appendText(
1027                                 "NOTE: The CPU has >25% more load than the GPU. This imbalance wastes\n"
1028                                 "      GPU resources.");
1029                     }
1030                 }
1031                 if (gpu_cpu_ratio > 1.25)
1032                 {
1033                     GMX_LOG(mdlog.warning).asParagraph().appendText(
1034                             "NOTE: The GPU has >25% more load than the CPU. This imbalance wastes\n"
1035                             "      CPU resources.");
1036                 }
1037             }
1038         }
1039     }
1040
1041     if (wc->wc_barrier)
1042     {
1043         GMX_LOG(mdlog.warning).asParagraph().appendText(
1044                 "MPI_Barrier was called before each cycle start/stop\n"
1045                 "call, so timings are not those of real runs.");
1046     }
1047
1048     if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
1049         (cyc_sum[ewcDOMDEC] > tot*0.1 ||
1050          cyc_sum[ewcNS] > tot*0.1))
1051     {
1052         /* Only the sim master calls this function, so always print to stderr */
1053         if (wc->wcc[ewcDOMDEC].n == 0)
1054         {
1055             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1056                     "NOTE: %d %% of the run time was spent in pair search,\n"
1057                     "      you might want to increase nstlist (this has no effect on accuracy)\n",
1058                     gmx::roundToInt(100*cyc_sum[ewcNS]/tot));
1059         }
1060         else
1061         {
1062             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1063                     "NOTE: %d %% of the run time was spent in domain decomposition,\n"
1064                     "      %d %% of the run time was spent in pair search,\n"
1065                     "      you might want to increase nstlist (this has no effect on accuracy)\n",
1066                     gmx::roundToInt(100*cyc_sum[ewcDOMDEC]/tot),
1067                     gmx::roundToInt(100*cyc_sum[ewcNS]/tot));
1068         }
1069     }
1070
1071     if (cyc_sum[ewcMoveE] > tot*0.05)
1072     {
1073         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1074                 "NOTE: %d %% of the run time was spent communicating energies,\n"
1075                 "      you might want to use the -gcom option of mdrun\n",
1076                 gmx::roundToInt(100*cyc_sum[ewcMoveE]/tot));
1077     }
1078 }
1079
1080 extern int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
1081 {
1082     if (wc == nullptr)
1083     {
1084         return -1;
1085     }
1086
1087     return wc->reset_counters;
1088 }
1089
1090 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, int64_t reset_counters)
1091 {
1092     if (wc == nullptr)
1093     {
1094         return;
1095     }
1096
1097     wc->reset_counters = reset_counters;
1098 }
1099
1100 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
1101 {
1102     if (useCycleSubcounters && wc != nullptr)
1103     {
1104         wc->wcsc[ewcs].start = gmx_cycles_read();
1105     }
1106 }
1107
1108 void wallcycle_sub_start_nocount(gmx_wallcycle_t wc, int ewcs)
1109 {
1110     if (useCycleSubcounters && wc != nullptr)
1111     {
1112         wallcycle_sub_start(wc, ewcs);
1113         wc->wcsc[ewcs].n--;
1114     }
1115 }
1116
1117 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
1118 {
1119     if (useCycleSubcounters && wc != nullptr)
1120     {
1121         wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
1122         wc->wcsc[ewcs].n++;
1123     }
1124 }