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