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