Remove unused function generate_excls and make clean_excls static
[gromacs.git] / src / gromacs / timing / wallcycle.cpp
blobcf2320591ecc470ddb8f54e0a83adeb00b925106
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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,2018,2019, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "wallcycle.h"
41 #include "config.h"
43 #include <cstdlib>
45 #include <array>
46 #include <vector>
48 #include "gromacs/math/functions.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/timing/cyclecounter.h"
51 #include "gromacs/timing/gpu_timing.h"
52 #include "gromacs/timing/wallcyclereporting.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/gmxmpi.h"
56 #include "gromacs/utility/logger.h"
57 #include "gromacs/utility/smalloc.h"
58 #include "gromacs/utility/snprintf.h"
60 static const bool useCycleSubcounters = GMX_CYCLE_SUBCOUNTERS;
62 /* DEBUG_WCYCLE adds consistency checking for the counters.
63 * It checks if you stop a counter different from the last
64 * one that was opened and if you do nest too deep.
66 /* #define DEBUG_WCYCLE */
68 #ifdef DEBUG_WCYCLE
69 #include "gromacs/utility/fatalerror.h"
70 #endif
72 typedef struct
74 int n;
75 gmx_cycles_t c;
76 gmx_cycles_t start;
77 } wallcc_t;
79 struct gmx_wallcycle
81 wallcc_t *wcc;
82 /* did we detect one or more invalid cycle counts */
83 gmx_bool haveInvalidCount;
84 /* variables for testing/debugging */
85 gmx_bool wc_barrier;
86 wallcc_t *wcc_all;
87 int wc_depth;
88 #ifdef DEBUG_WCYCLE
89 #define DEPTH_MAX 6
90 int counterlist[DEPTH_MAX];
91 int count_depth;
92 #endif
93 int ewc_prev;
94 gmx_cycles_t cycle_prev;
95 int64_t reset_counters;
96 #if GMX_MPI
97 MPI_Comm mpi_comm_mygroup;
98 #endif
99 wallcc_t *wcsc;
102 /* Each name should not exceed 19 printing characters
103 (ie. terminating null can be twentieth) */
104 static const char *wcn[ewcNR] =
106 "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
107 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
108 "Comm. coord.", "Force", "Wait + Comm. F", "PME mesh",
109 "PME redist. X/F", "PME spread", "PME gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve LJ", "PME solve Elec",
110 "PME wait for PP", "Wait + Recv. PME F",
111 "Wait PME GPU spread", "PME 3D-FFT", "PME solve", /* the strings for FFT/solve are repeated here for mixed mode counters */
112 "Wait PME GPU gather", "Wait Bonded GPU", "Reduce GPU PME F",
113 "Wait GPU NB nonloc.", "Wait GPU NB local", "NB X/F buffer ops.",
114 "Vsite spread", "COM pull force", "AWH",
115 "Write traj.", "Update", "Constraints", "Comm. energies",
116 "Enforced rotation", "Add rot. forces", "Position swapping", "IMD", "Test"
119 static const char *wcsn[ewcsNR] =
121 "DD redist.", "DD NS grid + sort", "DD setup comm.",
122 "DD make top.", "DD make constr.", "DD top. other",
123 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
124 "Bonded F",
125 "Bonded-FEP F",
126 "Restraints F",
127 "Listed buffer ops.",
128 "Nonbonded pruning",
129 "Nonbonded F kernel",
130 "Nonbonded F clear",
131 "Nonbonded FEP",
132 "Launch NB GPU tasks",
133 "Launch Bonded GPU tasks",
134 "Launch PME GPU tasks",
135 "Ewald F correction",
136 "NB X buffer ops.",
137 "NB F buffer ops.",
138 "Clear force buffer",
139 "Test subcounter",
142 /* PME GPU timing events' names - correspond to the enum in the gpu_timing.h */
143 static const char *PMEStageNames[] =
145 "PME spline",
146 "PME spread",
147 "PME spline + spread",
148 "PME 3D-FFT r2c",
149 "PME solve",
150 "PME 3D-FFT c2r",
151 "PME gather",
154 gmx_bool wallcycle_have_counter()
156 return gmx_cycles_have_counter();
159 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec gmx_unused *cr)
161 gmx_wallcycle_t wc;
164 if (!wallcycle_have_counter())
166 return nullptr;
169 snew(wc, 1);
171 wc->haveInvalidCount = FALSE;
172 wc->wc_barrier = FALSE;
173 wc->wcc_all = nullptr;
174 wc->wc_depth = 0;
175 wc->ewc_prev = -1;
176 wc->reset_counters = resetstep;
178 #if GMX_MPI
179 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != nullptr)
181 if (fplog)
183 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
185 wc->wc_barrier = TRUE;
186 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
188 #endif
190 snew(wc->wcc, ewcNR);
191 if (getenv("GMX_CYCLE_ALL") != nullptr)
193 if (fplog)
195 fprintf(fplog, "\nWill time all the code during the run\n\n");
197 snew(wc->wcc_all, ewcNR*ewcNR);
200 if (useCycleSubcounters)
202 snew(wc->wcsc, ewcsNR);
205 #ifdef DEBUG_WCYCLE
206 wc->count_depth = 0;
207 #endif
209 return wc;
212 void wallcycle_destroy(gmx_wallcycle_t wc)
214 if (wc == nullptr)
216 return;
219 if (wc->wcc != nullptr)
221 sfree(wc->wcc);
223 if (wc->wcc_all != nullptr)
225 sfree(wc->wcc_all);
227 if (wc->wcsc != nullptr)
229 sfree(wc->wcsc);
231 sfree(wc);
234 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
236 wc->ewc_prev = ewc;
237 wc->cycle_prev = cycle;
240 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
242 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
243 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
247 #ifdef DEBUG_WCYCLE
248 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
250 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
252 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
254 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
255 wc->count_depth);
257 wc->counterlist[wc->count_depth] = ewc;
258 wc->count_depth++;
261 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
263 wc->count_depth--;
265 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
267 if (wc->count_depth < 0)
269 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
271 if (wc->counterlist[wc->count_depth] != ewc)
273 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
274 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
277 #endif
279 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
281 gmx_cycles_t cycle;
283 if (wc == nullptr)
285 return;
288 #if GMX_MPI
289 if (wc->wc_barrier)
291 MPI_Barrier(wc->mpi_comm_mygroup);
293 #endif
295 #ifdef DEBUG_WCYCLE
296 debug_start_check(wc, ewc);
297 #endif
299 cycle = gmx_cycles_read();
300 wc->wcc[ewc].start = cycle;
301 if (wc->wcc_all != nullptr)
303 wc->wc_depth++;
304 if (ewc == ewcRUN)
306 wallcycle_all_start(wc, ewc, cycle);
308 else if (wc->wc_depth == 3)
310 wallcycle_all_stop(wc, ewc, cycle);
315 void wallcycle_increment_event_count(gmx_wallcycle_t wc, int ewc)
317 if (wc == nullptr)
319 return;
321 wc->wcc[ewc].n++;
324 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
326 if (wc == nullptr)
328 return;
331 wallcycle_start(wc, ewc);
332 wc->wcc[ewc].n--;
335 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
337 gmx_cycles_t cycle, last;
339 if (wc == nullptr)
341 return 0;
344 #if GMX_MPI
345 if (wc->wc_barrier)
347 MPI_Barrier(wc->mpi_comm_mygroup);
349 #endif
351 #ifdef DEBUG_WCYCLE
352 debug_stop_check(wc, ewc);
353 #endif
355 /* When processes or threads migrate between cores, the cycle counting
356 * can get messed up if the cycle counter on different cores are not
357 * synchronized. When this happens we expect both large negative and
358 * positive cycle differences. We can detect negative cycle differences.
359 * Detecting too large positive counts if difficult, since count can be
360 * large, especially for ewcRUN. If we detect a negative count,
361 * we will not print the cycle accounting table.
363 cycle = gmx_cycles_read();
364 if (cycle >= wc->wcc[ewc].start)
366 last = cycle - wc->wcc[ewc].start;
368 else
370 last = 0;
371 wc->haveInvalidCount = TRUE;
373 wc->wcc[ewc].c += last;
374 wc->wcc[ewc].n++;
375 if (wc->wcc_all)
377 wc->wc_depth--;
378 if (ewc == ewcRUN)
380 wallcycle_all_stop(wc, ewc, cycle);
382 else if (wc->wc_depth == 2)
384 wallcycle_all_start(wc, ewc, cycle);
388 return last;
391 void wallcycle_get(gmx_wallcycle_t wc, int ewc, int *n, double *c)
393 *n = wc->wcc[ewc].n;
394 *c = static_cast<double>(wc->wcc[ewc].c);
397 void wallcycle_reset_all(gmx_wallcycle_t wc)
399 int i;
401 if (wc == nullptr)
403 return;
406 for (i = 0; i < ewcNR; i++)
408 wc->wcc[i].n = 0;
409 wc->wcc[i].c = 0;
411 wc->haveInvalidCount = FALSE;
413 if (wc->wcc_all)
415 for (i = 0; i < ewcNR*ewcNR; i++)
417 wc->wcc_all[i].n = 0;
418 wc->wcc_all[i].c = 0;
421 if (wc->wcsc)
423 for (i = 0; i < ewcsNR; i++)
425 wc->wcsc[i].n = 0;
426 wc->wcsc[i].c = 0;
431 static gmx_bool is_pme_counter(int ewc)
433 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
436 static gmx_bool is_pme_subcounter(int ewc)
438 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
441 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
442 static void subtract_cycles(wallcc_t *wcc, int ewc_main, int ewc_sub)
444 if (wcc[ewc_sub].n > 0)
446 if (wcc[ewc_main].c >= wcc[ewc_sub].c)
448 wcc[ewc_main].c -= wcc[ewc_sub].c;
450 else
452 /* Something is wrong with the cycle counting */
453 wcc[ewc_main].c = 0;
458 void wallcycle_scale_by_num_threads(gmx_wallcycle_t wc, bool isPmeRank, int nthreads_pp, int nthreads_pme)
460 if (wc == nullptr)
462 return;
465 for (int i = 0; i < ewcNR; i++)
467 if (is_pme_counter(i) || (i == ewcRUN && isPmeRank))
469 wc->wcc[i].c *= nthreads_pme;
471 if (wc->wcc_all)
473 for (int j = 0; j < ewcNR; j++)
475 wc->wcc_all[i*ewcNR+j].c *= nthreads_pme;
479 else
481 wc->wcc[i].c *= nthreads_pp;
483 if (wc->wcc_all)
485 for (int j = 0; j < ewcNR; j++)
487 wc->wcc_all[i*ewcNR+j].c *= nthreads_pp;
492 if (useCycleSubcounters && wc->wcsc && !isPmeRank)
494 for (int i = 0; i < ewcsNR; i++)
496 wc->wcsc[i].c *= nthreads_pp;
501 /* TODO Make an object for this function to return, containing some
502 * vectors of something like wallcc_t for the summed wcc, wcc_all and
503 * wcsc, AND the original wcc for rank 0.
505 * The GPU timing is reported only for rank 0, so we want to preserve
506 * the original wcycle on that rank. Rank 0 also reports the global
507 * counts before that, so needs something to contain the global data
508 * without over-writing the rank-0 data. The current implementation
509 * uses cycles_sum to manage this, which works OK now because wcsc and
510 * wcc_all are unused by the GPU reporting, but it is not satisfactory
511 * for the future. Also, there's no need for MPI_Allreduce, since
512 * only MASTERRANK uses any of the results. */
513 WallcycleCounts wallcycle_sum(const t_commrec *cr, gmx_wallcycle_t wc)
515 WallcycleCounts cycles_sum;
516 wallcc_t *wcc;
517 double cycles[ewcNR+ewcsNR];
518 #if GMX_MPI
519 double cycles_n[ewcNR+ewcsNR+1];
520 #endif
521 int i;
522 int nsum;
524 if (wc == nullptr)
526 /* Default construction of std::array of non-class T can leave
527 the values indeterminate, just like a C array, and icc
528 warns about it. */
529 cycles_sum.fill(0);
530 return cycles_sum;
533 wcc = wc->wcc;
535 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
536 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
538 subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
540 if (cr->npmenodes == 0)
542 /* All nodes do PME (or no PME at all) */
543 subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
545 else
547 /* The are PME-only nodes */
548 if (wcc[ewcPMEMESH].n > 0)
550 /* This must be a PME only node, calculate the Wait + Comm. time */
551 GMX_ASSERT(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c, "Total run ticks must be greater than PME-only ticks");
552 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
556 /* Store the cycles in a double buffer for summing */
557 for (i = 0; i < ewcNR; i++)
559 #if GMX_MPI
560 cycles_n[i] = static_cast<double>(wcc[i].n);
561 #endif
562 cycles[i] = static_cast<double>(wcc[i].c);
564 nsum = ewcNR;
565 if (wc->wcsc)
567 for (i = 0; i < ewcsNR; i++)
569 #if GMX_MPI
570 cycles_n[ewcNR+i] = static_cast<double>(wc->wcsc[i].n);
571 #endif
572 cycles[ewcNR+i] = static_cast<double>(wc->wcsc[i].c);
574 nsum += ewcsNR;
577 #if GMX_MPI
578 if (cr->nnodes > 1)
580 double buf[ewcNR+ewcsNR+1];
582 // TODO this code is used only at the end of the run, so we
583 // can just do a simple reduce of haveInvalidCount in
584 // wallcycle_print, and avoid bugs
585 cycles_n[nsum] = (wc->haveInvalidCount ? 1 : 0);
586 // TODO Use MPI_Reduce
587 MPI_Allreduce(cycles_n, buf, nsum + 1, MPI_DOUBLE, MPI_MAX,
588 cr->mpi_comm_mysim);
589 for (i = 0; i < ewcNR; i++)
591 wcc[i].n = gmx::roundToInt(buf[i]);
593 wc->haveInvalidCount = (buf[nsum] > 0);
594 if (wc->wcsc)
596 for (i = 0; i < ewcsNR; i++)
598 wc->wcsc[i].n = gmx::roundToInt(buf[ewcNR+i]);
602 // TODO Use MPI_Reduce
603 MPI_Allreduce(cycles, cycles_sum.data(), nsum, MPI_DOUBLE, MPI_SUM,
604 cr->mpi_comm_mysim);
606 if (wc->wcc_all != nullptr)
608 double *buf_all, *cyc_all;
610 snew(cyc_all, ewcNR*ewcNR);
611 snew(buf_all, ewcNR*ewcNR);
612 for (i = 0; i < ewcNR*ewcNR; i++)
614 cyc_all[i] = wc->wcc_all[i].c;
616 // TODO Use MPI_Reduce
617 MPI_Allreduce(cyc_all, buf_all, ewcNR*ewcNR, MPI_DOUBLE, MPI_SUM,
618 cr->mpi_comm_mysim);
619 for (i = 0; i < ewcNR*ewcNR; i++)
621 wc->wcc_all[i].c = static_cast<gmx_cycles_t>(buf_all[i]);
623 sfree(buf_all);
624 sfree(cyc_all);
627 else
628 #endif
630 for (i = 0; i < nsum; i++)
632 cycles_sum[i] = cycles[i];
636 return cycles_sum;
639 static void print_cycles(FILE *fplog, double c2t, const char *name,
640 int nnodes, int nthreads,
641 int ncalls, double c_sum, double tot)
643 char nnodes_str[STRLEN];
644 char nthreads_str[STRLEN];
645 char ncalls_str[STRLEN];
646 double wallt;
647 double percentage = (tot > 0.) ? (100. * c_sum / tot) : 0.;
649 if (c_sum > 0)
651 if (ncalls > 0)
653 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
654 if (nnodes < 0)
656 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
658 else
660 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
662 if (nthreads < 0)
664 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
666 else
668 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
671 else
673 nnodes_str[0] = 0;
674 nthreads_str[0] = 0;
675 ncalls_str[0] = 0;
677 /* Convert the cycle count to wallclock time for this task */
678 wallt = c_sum*c2t;
680 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
681 name, nnodes_str, nthreads_str, ncalls_str, wallt,
682 c_sum*1e-9, percentage);
686 static void print_gputimes(FILE *fplog, const char *name,
687 int n, double t, double tot_t)
689 char num[11];
690 char avg_perf[11];
692 if (n > 0)
694 snprintf(num, sizeof(num), "%10d", n);
695 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
697 else
699 sprintf(num, " ");
700 sprintf(avg_perf, " ");
702 if (t != tot_t && tot_t > 0)
704 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
705 name, num, t/1000, avg_perf, 100 * t/tot_t);
707 else
709 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
710 name, "", t/1000, avg_perf, 100.0);
714 static void print_header(FILE *fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
716 int nrank_tot = nrank_pp + nrank_pme;
717 if (0 == nrank_pme)
719 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
720 if (nth_pp > 1)
722 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
724 /* Don't report doing PP+PME, because we can't tell here if
725 * this is RF, etc. */
727 else
729 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
730 if (nth_pp > 1)
732 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
734 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
735 if (nth_pme > 1)
737 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
741 fprintf(fplog, "\n\n");
742 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
743 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
747 void wallcycle_print(FILE *fplog, const gmx::MDLogger &mdlog, int nnodes, int npme,
748 int nth_pp, int nth_pme, double realtime,
749 gmx_wallcycle_t wc, const WallcycleCounts &cyc_sum,
750 const gmx_wallclock_gpu_nbnxn_t *gpu_nbnxn_t,
751 const gmx_wallclock_gpu_pme_t *gpu_pme_t)
753 double tot, tot_for_pp, tot_for_rest, tot_cpu_overlap, gpu_cpu_ratio;
754 double c2t, c2t_pp, c2t_pme = 0;
755 int i, j, npp, nth_tot;
756 char buf[STRLEN];
757 const char *hline = "-----------------------------------------------------------------------------";
759 if (wc == nullptr)
761 return;
764 GMX_ASSERT(nth_pp > 0, "Number of particle-particle threads must be >0");
765 GMX_ASSERT(nth_pme > 0, "Number of PME threads must be >0");
766 GMX_ASSERT(nnodes > 0, "Number of nodes must be >0");
767 GMX_ASSERT(npme >= 0, "Number of PME nodes cannot be negative");
768 npp = nnodes - npme;
769 /* npme is the number of PME-only ranks used, and we always do PP work */
770 GMX_ASSERT(npp > 0, "Number of particle-particle nodes must be >0");
772 nth_tot = npp*nth_pp + npme*nth_pme;
774 /* When using PME-only nodes, the next line is valid for both
775 PP-only and PME-only nodes because they started ewcRUN at the
776 same time. */
777 tot = cyc_sum[ewcRUN];
778 tot_for_pp = 0;
780 if (tot <= 0.0)
782 /* TODO This is heavy handed, but until someone reworks the
783 code so that it is provably robust with respect to
784 non-positive values for all possible timer and cycle
785 counters, there is less value gained from printing whatever
786 timing data might still be sensible for some non-Jenkins
787 run, than is lost from diagnosing Jenkins FP exceptions on
788 runs about whose execution time we don't care. */
789 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
790 "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a time accounting",
791 tot);
792 return;
795 if (wc->haveInvalidCount)
797 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.");
798 return;
802 /* Conversion factor from cycles to seconds */
803 c2t = realtime/tot;
804 c2t_pp = c2t * nth_tot / static_cast<double>(npp*nth_pp);
805 if (npme > 0)
807 c2t_pme = c2t * nth_tot / static_cast<double>(npme*nth_pme);
809 else
811 c2t_pme = 0;
814 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");
816 print_header(fplog, npp, nth_pp, npme, nth_pme);
818 fprintf(fplog, "%s\n", hline);
819 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
821 if (is_pme_subcounter(i))
823 /* Do not count these at all */
825 else if (npme > 0 && is_pme_counter(i))
827 /* Print timing information for PME-only nodes, but add an
828 * asterisk so the reader of the table can know that the
829 * walltimes are not meant to add up. The asterisk still
830 * fits in the required maximum of 19 characters. */
831 char buffer[STRLEN];
832 snprintf(buffer, STRLEN, "%s *", wcn[i]);
833 print_cycles(fplog, c2t_pme, buffer,
834 npme, nth_pme,
835 wc->wcc[i].n, cyc_sum[i], tot);
837 else
839 /* Print timing information when it is for a PP or PP+PME
840 node */
841 print_cycles(fplog, c2t_pp, wcn[i],
842 npp, nth_pp,
843 wc->wcc[i].n, cyc_sum[i], tot);
844 tot_for_pp += cyc_sum[i];
847 if (wc->wcc_all != nullptr)
849 for (i = 0; i < ewcNR; i++)
851 for (j = 0; j < ewcNR; j++)
853 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
854 print_cycles(fplog, c2t_pp, buf,
855 npp, nth_pp,
856 wc->wcc_all[i*ewcNR+j].n,
857 wc->wcc_all[i*ewcNR+j].c,
858 tot);
862 tot_for_rest = tot * npp * nth_pp / static_cast<double>(nth_tot);
863 print_cycles(fplog, c2t_pp, "Rest",
864 npp, nth_pp,
865 -1, tot_for_rest - tot_for_pp, tot);
866 fprintf(fplog, "%s\n", hline);
867 print_cycles(fplog, c2t, "Total",
868 npp, nth_pp,
869 -1, tot, tot);
870 fprintf(fplog, "%s\n", hline);
872 if (npme > 0)
874 fprintf(fplog,
875 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
876 " twice the total reported, but the cycle count total and %% are correct.\n"
877 "%s\n", hline);
880 if (wc->wcc[ewcPMEMESH].n > 0)
882 // A workaround to not print breakdown when no subcounters were recorded.
883 // TODO: figure out and record PME GPU counters (what to do with the waiting ones?)
884 std::vector<int> validPmeSubcounterIndices;
885 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
887 if (is_pme_subcounter(i) && wc->wcc[i].n > 0)
889 validPmeSubcounterIndices.push_back(i);
893 if (!validPmeSubcounterIndices.empty())
895 fprintf(fplog, " Breakdown of PME mesh computation\n");
896 fprintf(fplog, "%s\n", hline);
897 for (auto i : validPmeSubcounterIndices)
899 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i],
900 npme > 0 ? npme : npp, nth_pme,
901 wc->wcc[i].n, cyc_sum[i], tot);
903 fprintf(fplog, "%s\n", hline);
907 if (useCycleSubcounters && wc->wcsc)
909 fprintf(fplog, " Breakdown of PP computation\n");
910 fprintf(fplog, "%s\n", hline);
911 for (i = 0; i < ewcsNR; i++)
913 print_cycles(fplog, c2t_pp, wcsn[i],
914 npp, nth_pp,
915 wc->wcsc[i].n, cyc_sum[ewcNR+i], tot);
917 fprintf(fplog, "%s\n", hline);
920 /* print GPU timing summary */
921 double tot_gpu = 0.0;
922 if (gpu_pme_t)
924 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
926 tot_gpu += gpu_pme_t->timing[k].t;
929 if (gpu_nbnxn_t)
931 const char *k_log_str[2][2] = {
932 {"Nonbonded F kernel", "Nonbonded F+ene k."},
933 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
935 tot_gpu += gpu_nbnxn_t->pl_h2d_t + gpu_nbnxn_t->nb_h2d_t + gpu_nbnxn_t->nb_d2h_t;
937 /* add up the kernel timings */
938 for (i = 0; i < 2; i++)
940 for (j = 0; j < 2; j++)
942 tot_gpu += gpu_nbnxn_t->ktime[i][j].t;
945 tot_gpu += gpu_nbnxn_t->pruneTime.t;
947 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
948 if (wc->wcc[ewcPMEMESH].n > 0)
950 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
952 tot_cpu_overlap *= realtime*1000/tot; /* convert s to ms */
954 fprintf(fplog, "\n GPU timings\n%s\n", hline);
955 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
956 fprintf(fplog, "%s\n", hline);
957 print_gputimes(fplog, "Pair list H2D",
958 gpu_nbnxn_t->pl_h2d_c, gpu_nbnxn_t->pl_h2d_t, tot_gpu);
959 print_gputimes(fplog, "X / q H2D",
960 gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_h2d_t, tot_gpu);
962 for (i = 0; i < 2; i++)
964 for (j = 0; j < 2; j++)
966 if (gpu_nbnxn_t->ktime[i][j].c)
968 print_gputimes(fplog, k_log_str[i][j],
969 gpu_nbnxn_t->ktime[i][j].c, gpu_nbnxn_t->ktime[i][j].t, tot_gpu);
973 if (gpu_pme_t)
975 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
977 if (gpu_pme_t->timing[k].c)
979 print_gputimes(fplog, PMEStageNames[k],
980 gpu_pme_t->timing[k].c,
981 gpu_pme_t->timing[k].t,
982 tot_gpu);
986 if (gpu_nbnxn_t->pruneTime.c)
988 print_gputimes(fplog, "Pruning kernel", gpu_nbnxn_t->pruneTime.c, gpu_nbnxn_t->pruneTime.t, tot_gpu);
990 print_gputimes(fplog, "F D2H", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_d2h_t, tot_gpu);
991 fprintf(fplog, "%s\n", hline);
992 print_gputimes(fplog, "Total ", gpu_nbnxn_t->nb_c, tot_gpu, tot_gpu);
993 fprintf(fplog, "%s\n", hline);
994 if (gpu_nbnxn_t->dynamicPruneTime.c)
996 /* We print the dynamic pruning kernel timings after a separator
997 * and avoid adding it to tot_gpu as this is not in the force
998 * overlap. We print the fraction as relative to the rest.
1000 print_gputimes(fplog, "*Dynamic pruning", gpu_nbnxn_t->dynamicPruneTime.c, gpu_nbnxn_t->dynamicPruneTime.t, tot_gpu);
1001 fprintf(fplog, "%s\n", hline);
1003 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
1004 if (gpu_nbnxn_t->nb_c > 0 && wc->wcc[ewcFORCE].n > 0)
1006 fprintf(fplog, "\nAverage per-step force GPU/CPU evaluation time ratio: %.3f ms/%.3f ms = %.3f\n",
1007 tot_gpu/gpu_nbnxn_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
1008 gpu_cpu_ratio);
1011 /* only print notes related to CPU-GPU load balance with PME */
1012 if (wc->wcc[ewcPMEMESH].n > 0)
1014 fprintf(fplog, "For optimal resource utilization this ratio should be close to 1\n");
1016 /* print note if the imbalance is high with PME case in which
1017 * CPU-GPU load balancing is possible */
1018 if (gpu_cpu_ratio < 0.8 || gpu_cpu_ratio > 1.25)
1020 /* Only the sim master calls this function, so always print to stderr */
1021 if (gpu_cpu_ratio < 0.8)
1023 if (npp > 1)
1025 /* The user could have used -notunepme,
1026 * but we currently can't check that here.
1028 GMX_LOG(mdlog.warning).asParagraph().appendText(
1029 "NOTE: The CPU has >25% more load than the GPU. This imbalance wastes\n"
1030 " GPU resources. Maybe the domain decomposition limits the PME tuning.\n"
1031 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
1033 else
1035 /* We should not end up here, unless the box is
1036 * too small for increasing the cut-off for PME tuning.
1038 GMX_LOG(mdlog.warning).asParagraph().appendText(
1039 "NOTE: The CPU has >25% more load than the GPU. This imbalance wastes\n"
1040 " GPU resources.");
1043 if (gpu_cpu_ratio > 1.25)
1045 GMX_LOG(mdlog.warning).asParagraph().appendText(
1046 "NOTE: The GPU has >25% more load than the CPU. This imbalance wastes\n"
1047 " CPU resources.");
1053 if (wc->wc_barrier)
1055 GMX_LOG(mdlog.warning).asParagraph().appendText(
1056 "MPI_Barrier was called before each cycle start/stop\n"
1057 "call, so timings are not those of real runs.");
1060 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
1061 (cyc_sum[ewcDOMDEC] > tot*0.1 ||
1062 cyc_sum[ewcNS] > tot*0.1))
1064 /* Only the sim master calls this function, so always print to stderr */
1065 if (wc->wcc[ewcDOMDEC].n == 0)
1067 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1068 "NOTE: %d %% of the run time was spent in pair search,\n"
1069 " you might want to increase nstlist (this has no effect on accuracy)\n",
1070 gmx::roundToInt(100*cyc_sum[ewcNS]/tot));
1072 else
1074 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1075 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
1076 " %d %% of the run time was spent in pair search,\n"
1077 " you might want to increase nstlist (this has no effect on accuracy)\n",
1078 gmx::roundToInt(100*cyc_sum[ewcDOMDEC]/tot),
1079 gmx::roundToInt(100*cyc_sum[ewcNS]/tot));
1083 if (cyc_sum[ewcMoveE] > tot*0.05)
1085 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1086 "NOTE: %d %% of the run time was spent communicating energies,\n"
1087 " you might want to increase some nst* mdp options\n",
1088 gmx::roundToInt(100*cyc_sum[ewcMoveE]/tot));
1092 extern int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
1094 if (wc == nullptr)
1096 return -1;
1099 return wc->reset_counters;
1102 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, int64_t reset_counters)
1104 if (wc == nullptr)
1106 return;
1109 wc->reset_counters = reset_counters;
1112 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
1114 if (useCycleSubcounters && wc != nullptr)
1116 wc->wcsc[ewcs].start = gmx_cycles_read();
1120 void wallcycle_sub_start_nocount(gmx_wallcycle_t wc, int ewcs)
1122 if (useCycleSubcounters && wc != nullptr)
1124 wallcycle_sub_start(wc, ewcs);
1125 wc->wcsc[ewcs].n--;
1129 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
1131 if (useCycleSubcounters && wc != nullptr)
1133 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
1134 wc->wcsc[ewcs].n++;