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