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