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