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