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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 #include "wallcycle.h"
49 #include "gromacs/math/functions.h"
50 #include "gromacs/mdtypes/commrec.h"
51 #include "gromacs/timing/cyclecounter.h"
52 #include "gromacs/timing/gpu_timing.h"
53 #include "gromacs/timing/wallcyclereporting.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/gmxassert.h"
56 #include "gromacs/utility/gmxmpi.h"
57 #include "gromacs/utility/logger.h"
58 #include "gromacs/utility/smalloc.h"
59 #include "gromacs/utility/snprintf.h"
61 static const bool useCycleSubcounters = GMX_CYCLE_SUBCOUNTERS;
63 /* DEBUG_WCYCLE adds consistency checking for the counters.
64 * It checks if you stop a counter different from the last
65 * one that was opened and if you do nest too deep.
67 /* #define DEBUG_WCYCLE */
70 # include "gromacs/utility/fatalerror.h"
83 /* did we detect one or more invalid cycle counts */
84 gmx_bool haveInvalidCount;
85 /* variables for testing/debugging */
91 int counterlist[DEPTH_MAX];
95 gmx_cycles_t cycle_prev;
96 int64_t reset_counters;
98 MPI_Comm mpi_comm_mygroup;
103 /* Each name should not exceed 19 printing characters
104 (ie. terminating null can be twentieth) */
105 static const char* wcn[ewcNR] = { "Run",
127 "Wait + Recv. PME F",
128 "Wait PME GPU spread",
130 "PME solve", /* the strings for FFT/solve are repeated here for mixed mode counters */
131 "Wait PME GPU gather",
134 "Wait GPU NB nonloc.",
136 "Wait GPU state copy",
137 "NB X/F buffer ops.",
151 static const char* wcsn[ewcsNR] = {
161 "NS search non-loc.",
165 "Listed buffer ops.",
167 "Nonbonded F kernel",
170 "Launch NB GPU tasks",
171 "Launch Bonded GPU tasks",
172 "Launch PME GPU tasks",
174 "Ewald F correction",
177 "Clear force buffer",
181 /* PME GPU timing events' names - correspond to the enum in the gpu_timing.h */
182 static const char* PMEStageNames[] = {
183 "PME spline", "PME spread", "PME spline + spread", "PME 3D-FFT r2c",
184 "PME solve", "PME 3D-FFT c2r", "PME gather",
187 gmx_bool wallcycle_have_counter()
189 return gmx_cycles_have_counter();
192 gmx_wallcycle_t wallcycle_init(FILE* fplog, int resetstep, t_commrec gmx_unused* cr)
197 if (!wallcycle_have_counter())
204 wc->haveInvalidCount = FALSE;
205 wc->wc_barrier = FALSE;
206 wc->wcc_all = nullptr;
209 wc->reset_counters = resetstep;
212 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != nullptr)
216 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
218 wc->wc_barrier = TRUE;
219 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
223 snew(wc->wcc, ewcNR);
224 if (getenv("GMX_CYCLE_ALL") != nullptr)
228 fprintf(fplog, "\nWill time all the code during the run\n\n");
230 snew(wc->wcc_all, ewcNR * ewcNR);
233 if (useCycleSubcounters)
235 snew(wc->wcsc, ewcsNR);
245 void wallcycle_destroy(gmx_wallcycle_t wc)
252 if (wc->wcc != nullptr)
256 if (wc->wcc_all != nullptr)
260 if (wc->wcsc != nullptr)
267 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
270 wc->cycle_prev = cycle;
273 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
275 wc->wcc_all[wc->ewc_prev * ewcNR + ewc].n += 1;
276 wc->wcc_all[wc->ewc_prev * ewcNR + ewc].c += cycle - wc->cycle_prev;
281 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
283 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
285 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
287 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d", wc->count_depth);
289 wc->counterlist[wc->count_depth] = ewc;
293 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
297 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
299 if (wc->count_depth < 0)
301 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc],
304 if (wc->counterlist[wc->count_depth] != ewc)
306 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
307 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
312 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
324 MPI_Barrier(wc->mpi_comm_mygroup);
329 debug_start_check(wc, ewc);
332 cycle = gmx_cycles_read();
333 wc->wcc[ewc].start = cycle;
334 if (wc->wcc_all != nullptr)
339 wallcycle_all_start(wc, ewc, cycle);
341 else if (wc->wc_depth == 3)
343 wallcycle_all_stop(wc, ewc, cycle);
348 void wallcycle_increment_event_count(gmx_wallcycle_t wc, int ewc)
357 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
364 wallcycle_start(wc, ewc);
368 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
370 gmx_cycles_t cycle, last;
380 MPI_Barrier(wc->mpi_comm_mygroup);
385 debug_stop_check(wc, ewc);
388 /* When processes or threads migrate between cores, the cycle counting
389 * can get messed up if the cycle counter on different cores are not
390 * synchronized. When this happens we expect both large negative and
391 * positive cycle differences. We can detect negative cycle differences.
392 * Detecting too large positive counts if difficult, since count can be
393 * large, especially for ewcRUN. If we detect a negative count,
394 * we will not print the cycle accounting table.
396 cycle = gmx_cycles_read();
397 if (cycle >= wc->wcc[ewc].start)
399 last = cycle - wc->wcc[ewc].start;
404 wc->haveInvalidCount = TRUE;
406 wc->wcc[ewc].c += last;
413 wallcycle_all_stop(wc, ewc, cycle);
415 else if (wc->wc_depth == 2)
417 wallcycle_all_start(wc, ewc, cycle);
424 void wallcycle_get(gmx_wallcycle_t wc, int ewc, int* n, double* c)
427 *c = static_cast<double>(wc->wcc[ewc].c);
430 void wallcycle_reset_all(gmx_wallcycle_t wc)
439 for (i = 0; i < ewcNR; i++)
444 wc->haveInvalidCount = FALSE;
448 for (i = 0; i < ewcNR * ewcNR; i++)
450 wc->wcc_all[i].n = 0;
451 wc->wcc_all[i].c = 0;
456 for (i = 0; i < ewcsNR; i++)
464 static gmx_bool is_pme_counter(int ewc)
466 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
469 static gmx_bool is_pme_subcounter(int ewc)
471 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
474 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
475 static void subtract_cycles(wallcc_t* wcc, int ewc_main, int ewc_sub)
477 if (wcc[ewc_sub].n > 0)
479 if (wcc[ewc_main].c >= wcc[ewc_sub].c)
481 wcc[ewc_main].c -= wcc[ewc_sub].c;
485 /* Something is wrong with the cycle counting */
491 void wallcycle_scale_by_num_threads(gmx_wallcycle_t wc, bool isPmeRank, int nthreads_pp, int nthreads_pme)
498 for (int i = 0; i < ewcNR; i++)
500 if (is_pme_counter(i) || (i == ewcRUN && isPmeRank))
502 wc->wcc[i].c *= nthreads_pme;
506 for (int j = 0; j < ewcNR; j++)
508 wc->wcc_all[i * ewcNR + j].c *= nthreads_pme;
514 wc->wcc[i].c *= nthreads_pp;
518 for (int j = 0; j < ewcNR; j++)
520 wc->wcc_all[i * ewcNR + j].c *= nthreads_pp;
525 if (useCycleSubcounters && wc->wcsc && !isPmeRank)
527 for (int i = 0; i < ewcsNR; i++)
529 wc->wcsc[i].c *= nthreads_pp;
534 /* TODO Make an object for this function to return, containing some
535 * vectors of something like wallcc_t for the summed wcc, wcc_all and
536 * wcsc, AND the original wcc for rank 0.
538 * The GPU timing is reported only for rank 0, so we want to preserve
539 * the original wcycle on that rank. Rank 0 also reports the global
540 * counts before that, so needs something to contain the global data
541 * without over-writing the rank-0 data. The current implementation
542 * uses cycles_sum to manage this, which works OK now because wcsc and
543 * wcc_all are unused by the GPU reporting, but it is not satisfactory
544 * for the future. Also, there's no need for MPI_Allreduce, since
545 * only MASTERRANK uses any of the results. */
546 WallcycleCounts wallcycle_sum(const t_commrec* cr, gmx_wallcycle_t wc)
548 WallcycleCounts cycles_sum;
550 double cycles[int(ewcNR) + int(ewcsNR)];
552 double cycles_n[int(ewcNR) + int(ewcsNR) + 1];
559 /* Default construction of std::array of non-class T can leave
560 the values indeterminate, just like a C array, and icc
568 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
569 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
571 subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
573 if (cr->npmenodes == 0)
575 /* All nodes do PME (or no PME at all) */
576 subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
580 /* The are PME-only nodes */
581 if (wcc[ewcPMEMESH].n > 0)
583 /* This must be a PME only node, calculate the Wait + Comm. time */
584 GMX_ASSERT(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c,
585 "Total run ticks must be greater than PME-only ticks");
586 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
590 /* Store the cycles in a double buffer for summing */
591 for (i = 0; i < ewcNR; i++)
594 cycles_n[i] = static_cast<double>(wcc[i].n);
596 cycles[i] = static_cast<double>(wcc[i].c);
601 for (i = 0; i < ewcsNR; i++)
604 cycles_n[ewcNR + i] = static_cast<double>(wc->wcsc[i].n);
606 cycles[ewcNR + i] = static_cast<double>(wc->wcsc[i].c);
614 double buf[int(ewcNR) + int(ewcsNR) + 1];
616 // TODO this code is used only at the end of the run, so we
617 // can just do a simple reduce of haveInvalidCount in
618 // wallcycle_print, and avoid bugs
619 cycles_n[nsum] = (wc->haveInvalidCount ? 1 : 0);
620 // TODO Use MPI_Reduce
621 MPI_Allreduce(cycles_n, buf, nsum + 1, MPI_DOUBLE, MPI_MAX, cr->mpi_comm_mysim);
622 for (i = 0; i < ewcNR; i++)
624 wcc[i].n = gmx::roundToInt(buf[i]);
626 wc->haveInvalidCount = (buf[nsum] > 0);
629 for (i = 0; i < ewcsNR; i++)
631 wc->wcsc[i].n = gmx::roundToInt(buf[ewcNR + i]);
635 // TODO Use MPI_Reduce
636 MPI_Allreduce(cycles, cycles_sum.data(), nsum, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
638 if (wc->wcc_all != nullptr)
640 double *buf_all, *cyc_all;
642 snew(cyc_all, ewcNR * ewcNR);
643 snew(buf_all, ewcNR * ewcNR);
644 for (i = 0; i < ewcNR * ewcNR; i++)
646 cyc_all[i] = wc->wcc_all[i].c;
648 // TODO Use MPI_Reduce
649 MPI_Allreduce(cyc_all, buf_all, ewcNR * ewcNR, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
650 for (i = 0; i < ewcNR * ewcNR; i++)
652 wc->wcc_all[i].c = static_cast<gmx_cycles_t>(buf_all[i]);
661 for (i = 0; i < nsum; i++)
663 cycles_sum[i] = cycles[i];
671 print_cycles(FILE* fplog, double c2t, const char* name, int nnodes, int nthreads, int ncalls, double c_sum, double tot)
673 char nnodes_str[STRLEN];
674 char nthreads_str[STRLEN];
675 char ncalls_str[STRLEN];
677 double percentage = (tot > 0.) ? (100. * c_sum / tot) : 0.;
683 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
686 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
690 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
694 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
698 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
707 /* Convert the cycle count to wallclock time for this task */
710 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n", name, nnodes_str,
711 nthreads_str, ncalls_str, wallt, c_sum * 1e-9, percentage);
715 static void print_gputimes(FILE* fplog, const char* name, int n, double t, double tot_t)
722 snprintf(num, sizeof(num), "%10d", n);
723 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t / n);
728 sprintf(avg_perf, " ");
730 if (t != tot_t && tot_t > 0)
732 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n", name, num, t / 1000, avg_perf, 100 * t / tot_t);
736 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n", name, "", t / 1000, avg_perf, 100.0);
740 static void print_header(FILE* fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
742 int nrank_tot = nrank_pp + nrank_pme;
745 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
748 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
750 /* Don't report doing PP+PME, because we can't tell here if
751 * this is RF, etc. */
755 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
758 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
760 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
763 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
767 fprintf(fplog, "\n\n");
768 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
769 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
773 void wallcycle_print(FILE* fplog,
774 const gmx::MDLogger& mdlog,
781 const WallcycleCounts& cyc_sum,
782 const gmx_wallclock_gpu_nbnxn_t* gpu_nbnxn_t,
783 const gmx_wallclock_gpu_pme_t* gpu_pme_t)
785 double tot, tot_for_pp, tot_for_rest, tot_cpu_overlap, gpu_cpu_ratio;
786 double c2t, c2t_pp, c2t_pme = 0;
787 int i, j, npp, nth_tot;
790 "-----------------------------------------------------------------------------";
797 GMX_ASSERT(nth_pp > 0, "Number of particle-particle threads must be >0");
798 GMX_ASSERT(nth_pme > 0, "Number of PME threads must be >0");
799 GMX_ASSERT(nnodes > 0, "Number of nodes must be >0");
800 GMX_ASSERT(npme >= 0, "Number of PME nodes cannot be negative");
802 /* npme is the number of PME-only ranks used, and we always do PP work */
803 GMX_ASSERT(npp > 0, "Number of particle-particle nodes must be >0");
805 nth_tot = npp * nth_pp + npme * nth_pme;
807 /* When using PME-only nodes, the next line is valid for both
808 PP-only and PME-only nodes because they started ewcRUN at the
810 tot = cyc_sum[ewcRUN];
815 /* TODO This is heavy handed, but until someone reworks the
816 code so that it is provably robust with respect to
817 non-positive values for all possible timer and cycle
818 counters, there is less value gained from printing whatever
819 timing data might still be sensible for some non-Jenkins
820 run, than is lost from diagnosing Jenkins FP exceptions on
821 runs about whose execution time we don't care. */
822 GMX_LOG(mdlog.warning)
824 .appendTextFormatted(
825 "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a "
831 if (wc->haveInvalidCount)
833 GMX_LOG(mdlog.warning)
836 "NOTE: Detected invalid cycle counts, probably because threads moved "
837 "between CPU cores that do not have synchronized cycle counters. Will not "
838 "print the cycle accounting.");
843 /* Conversion factor from cycles to seconds */
844 c2t = realtime / tot;
845 c2t_pp = c2t * nth_tot / static_cast<double>(npp * nth_pp);
848 c2t_pme = c2t * nth_tot / static_cast<double>(npme * nth_pme);
855 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");
857 print_header(fplog, npp, nth_pp, npme, nth_pme);
859 fprintf(fplog, "%s\n", hline);
860 for (i = ewcPPDURINGPME + 1; i < ewcNR; i++)
862 if (is_pme_subcounter(i))
864 /* Do not count these at all */
866 else if (npme > 0 && is_pme_counter(i))
868 /* Print timing information for PME-only nodes, but add an
869 * asterisk so the reader of the table can know that the
870 * walltimes are not meant to add up. The asterisk still
871 * fits in the required maximum of 19 characters. */
873 snprintf(buffer, STRLEN, "%s *", wcn[i]);
874 print_cycles(fplog, c2t_pme, buffer, npme, nth_pme, wc->wcc[i].n, cyc_sum[i], tot);
878 /* Print timing information when it is for a PP or PP+PME
880 print_cycles(fplog, c2t_pp, wcn[i], npp, nth_pp, wc->wcc[i].n, cyc_sum[i], tot);
881 tot_for_pp += cyc_sum[i];
884 if (wc->wcc_all != nullptr)
886 for (i = 0; i < ewcNR; i++)
888 for (j = 0; j < ewcNR; j++)
890 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
891 print_cycles(fplog, c2t_pp, buf, npp, nth_pp, wc->wcc_all[i * ewcNR + j].n,
892 wc->wcc_all[i * ewcNR + j].c, tot);
896 tot_for_rest = tot * npp * nth_pp / static_cast<double>(nth_tot);
897 print_cycles(fplog, c2t_pp, "Rest", npp, nth_pp, -1, tot_for_rest - tot_for_pp, tot);
898 fprintf(fplog, "%s\n", hline);
899 print_cycles(fplog, c2t, "Total", npp, nth_pp, -1, tot, tot);
900 fprintf(fplog, "%s\n", hline);
905 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
906 " twice the total reported, but the cycle count total and %% are correct.\n"
911 if (wc->wcc[ewcPMEMESH].n > 0)
913 // A workaround to not print breakdown when no subcounters were recorded.
914 // TODO: figure out and record PME GPU counters (what to do with the waiting ones?)
915 std::vector<int> validPmeSubcounterIndices;
916 for (i = ewcPPDURINGPME + 1; i < ewcNR; i++)
918 if (is_pme_subcounter(i) && wc->wcc[i].n > 0)
920 validPmeSubcounterIndices.push_back(i);
924 if (!validPmeSubcounterIndices.empty())
926 fprintf(fplog, " Breakdown of PME mesh computation\n");
927 fprintf(fplog, "%s\n", hline);
928 for (auto i : validPmeSubcounterIndices)
930 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i], npme > 0 ? npme : npp,
931 nth_pme, wc->wcc[i].n, cyc_sum[i], tot);
933 fprintf(fplog, "%s\n", hline);
937 if (useCycleSubcounters && wc->wcsc)
939 fprintf(fplog, " Breakdown of PP computation\n");
940 fprintf(fplog, "%s\n", hline);
941 for (i = 0; i < ewcsNR; i++)
943 print_cycles(fplog, c2t_pp, wcsn[i], npp, nth_pp, wc->wcsc[i].n, cyc_sum[ewcNR + i], tot);
945 fprintf(fplog, "%s\n", hline);
948 /* print GPU timing summary */
949 double tot_gpu = 0.0;
952 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
954 tot_gpu += gpu_pme_t->timing[k].t;
959 const char* k_log_str[2][2] = { { "Nonbonded F kernel", "Nonbonded F+ene k." },
960 { "Nonbonded F+prune k.", "Nonbonded F+ene+prune k." } };
961 tot_gpu += gpu_nbnxn_t->pl_h2d_t + gpu_nbnxn_t->nb_h2d_t + gpu_nbnxn_t->nb_d2h_t;
963 /* add up the kernel timings */
964 for (i = 0; i < 2; i++)
966 for (j = 0; j < 2; j++)
968 tot_gpu += gpu_nbnxn_t->ktime[i][j].t;
971 tot_gpu += gpu_nbnxn_t->pruneTime.t;
973 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
974 if (wc->wcc[ewcPMEMESH].n > 0)
976 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
978 tot_cpu_overlap *= realtime * 1000 / tot; /* convert s to ms */
980 fprintf(fplog, "\n GPU timings\n%s\n", hline);
982 " Computing: Count Wall t (s) ms/step %c\n", '%');
983 fprintf(fplog, "%s\n", hline);
984 print_gputimes(fplog, "Pair list H2D", gpu_nbnxn_t->pl_h2d_c, gpu_nbnxn_t->pl_h2d_t, tot_gpu);
985 print_gputimes(fplog, "X / q H2D", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_h2d_t, tot_gpu);
987 for (i = 0; i < 2; i++)
989 for (j = 0; j < 2; j++)
991 if (gpu_nbnxn_t->ktime[i][j].c)
993 print_gputimes(fplog, k_log_str[i][j], gpu_nbnxn_t->ktime[i][j].c,
994 gpu_nbnxn_t->ktime[i][j].t, tot_gpu);
1000 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
1002 if (gpu_pme_t->timing[k].c)
1004 print_gputimes(fplog, PMEStageNames[k], gpu_pme_t->timing[k].c,
1005 gpu_pme_t->timing[k].t, tot_gpu);
1009 if (gpu_nbnxn_t->pruneTime.c)
1011 print_gputimes(fplog, "Pruning kernel", gpu_nbnxn_t->pruneTime.c,
1012 gpu_nbnxn_t->pruneTime.t, tot_gpu);
1014 print_gputimes(fplog, "F D2H", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_d2h_t, tot_gpu);
1015 fprintf(fplog, "%s\n", hline);
1016 print_gputimes(fplog, "Total ", gpu_nbnxn_t->nb_c, tot_gpu, tot_gpu);
1017 fprintf(fplog, "%s\n", hline);
1018 if (gpu_nbnxn_t->dynamicPruneTime.c)
1020 /* We print the dynamic pruning kernel timings after a separator
1021 * and avoid adding it to tot_gpu as this is not in the force
1022 * overlap. We print the fraction as relative to the rest.
1024 print_gputimes(fplog, "*Dynamic pruning", gpu_nbnxn_t->dynamicPruneTime.c,
1025 gpu_nbnxn_t->dynamicPruneTime.t, tot_gpu);
1026 fprintf(fplog, "%s\n", hline);
1028 gpu_cpu_ratio = tot_gpu / tot_cpu_overlap;
1029 if (gpu_nbnxn_t->nb_c > 0 && wc->wcc[ewcFORCE].n > 0)
1032 "\nAverage per-step force GPU/CPU evaluation time ratio: %.3f ms/%.3f ms = "
1034 tot_gpu / gpu_nbnxn_t->nb_c, tot_cpu_overlap / wc->wcc[ewcFORCE].n, gpu_cpu_ratio);
1037 /* only print notes related to CPU-GPU load balance with PME */
1038 if (wc->wcc[ewcPMEMESH].n > 0)
1040 fprintf(fplog, "For optimal resource utilization this ratio should be close to 1\n");
1042 /* print note if the imbalance is high with PME case in which
1043 * CPU-GPU load balancing is possible */
1044 if (gpu_cpu_ratio < 0.8 || gpu_cpu_ratio > 1.25)
1046 /* Only the sim master calls this function, so always print to stderr */
1047 if (gpu_cpu_ratio < 0.8)
1051 /* The user could have used -notunepme,
1052 * but we currently can't check that here.
1054 GMX_LOG(mdlog.warning)
1057 "NOTE: The CPU has >25% more load than the GPU. This "
1058 "imbalance wastes\n"
1059 " GPU resources. Maybe the domain decomposition "
1060 "limits the PME tuning.\n"
1061 " In that case, try setting the DD grid manually "
1062 "(-dd) or lowering -dds.");
1066 /* We should not end up here, unless the box is
1067 * too small for increasing the cut-off for PME tuning.
1069 GMX_LOG(mdlog.warning)
1072 "NOTE: The CPU has >25% more load than the GPU. This "
1073 "imbalance wastes\n"
1077 if (gpu_cpu_ratio > 1.25)
1079 GMX_LOG(mdlog.warning)
1082 "NOTE: The GPU has >25% more load than the CPU. This imbalance "
1092 GMX_LOG(mdlog.warning)
1095 "MPI_Barrier was called before each cycle start/stop\n"
1096 "call, so timings are not those of real runs.");
1099 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 && (cyc_sum[ewcDOMDEC] > tot * 0.1 || cyc_sum[ewcNS] > tot * 0.1))
1101 /* Only the sim master calls this function, so always print to stderr */
1102 if (wc->wcc[ewcDOMDEC].n == 0)
1104 GMX_LOG(mdlog.warning)
1106 .appendTextFormatted(
1107 "NOTE: %d %% of the run time was spent in pair search,\n"
1108 " you might want to increase nstlist (this has no effect on "
1110 gmx::roundToInt(100 * cyc_sum[ewcNS] / tot));
1114 GMX_LOG(mdlog.warning)
1116 .appendTextFormatted(
1117 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
1118 " %d %% of the run time was spent in pair search,\n"
1119 " you might want to increase nstlist (this has no effect on "
1121 gmx::roundToInt(100 * cyc_sum[ewcDOMDEC] / tot),
1122 gmx::roundToInt(100 * cyc_sum[ewcNS] / tot));
1126 if (cyc_sum[ewcMoveE] > tot * 0.05)
1128 GMX_LOG(mdlog.warning)
1130 .appendTextFormatted(
1131 "NOTE: %d %% of the run time was spent communicating energies,\n"
1132 " you might want to increase some nst* mdp options\n",
1133 gmx::roundToInt(100 * cyc_sum[ewcMoveE] / tot));
1137 extern int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
1144 return wc->reset_counters;
1147 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, int64_t reset_counters)
1154 wc->reset_counters = reset_counters;
1157 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
1159 if (useCycleSubcounters && wc != nullptr)
1161 wc->wcsc[ewcs].start = gmx_cycles_read();
1165 void wallcycle_sub_start_nocount(gmx_wallcycle_t wc, int ewcs)
1167 if (useCycleSubcounters && wc != nullptr)
1169 wallcycle_sub_start(wc, ewcs);
1174 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
1176 if (useCycleSubcounters && wc != nullptr)
1178 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;