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