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