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