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.
39 #include "wallcycle.h"
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 */
69 # include "gromacs/utility/fatalerror.h"
82 /* did we detect one or more invalid cycle counts */
83 gmx_bool haveInvalidCount;
84 /* variables for testing/debugging */
90 int counterlist[DEPTH_MAX];
94 gmx_cycles_t cycle_prev;
95 int64_t reset_counters;
97 MPI_Comm mpi_comm_mygroup;
102 /* Each name should not exceed 19 printing characters
103 (ie. terminating null can be twentieth) */
104 static const char* wcn[ewcNR] = { "Run",
126 "Wait + Recv. PME F",
127 "Wait PME GPU spread",
129 "PME solve", /* the strings for FFT/solve are repeated here for mixed mode counters */
130 "Wait PME GPU gather",
133 "Wait GPU NB nonloc.",
135 "NB X/F buffer ops.",
149 static const char* wcsn[ewcsNR] = {
159 "NS search non-loc.",
163 "Listed buffer ops.",
165 "Nonbonded F kernel",
168 "Launch NB GPU tasks",
169 "Launch Bonded GPU tasks",
170 "Launch PME GPU tasks",
171 "Ewald F correction",
174 "Clear force buffer",
178 /* PME GPU timing events' names - correspond to the enum in the gpu_timing.h */
179 static const char* PMEStageNames[] = {
180 "PME spline", "PME spread", "PME spline + spread", "PME 3D-FFT r2c",
181 "PME solve", "PME 3D-FFT c2r", "PME gather",
184 gmx_bool wallcycle_have_counter()
186 return gmx_cycles_have_counter();
189 gmx_wallcycle_t wallcycle_init(FILE* fplog, int resetstep, t_commrec gmx_unused* cr)
194 if (!wallcycle_have_counter())
201 wc->haveInvalidCount = FALSE;
202 wc->wc_barrier = FALSE;
203 wc->wcc_all = nullptr;
206 wc->reset_counters = resetstep;
209 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != nullptr)
213 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
215 wc->wc_barrier = TRUE;
216 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
220 snew(wc->wcc, ewcNR);
221 if (getenv("GMX_CYCLE_ALL") != nullptr)
225 fprintf(fplog, "\nWill time all the code during the run\n\n");
227 snew(wc->wcc_all, ewcNR * ewcNR);
230 if (useCycleSubcounters)
232 snew(wc->wcsc, ewcsNR);
242 void wallcycle_destroy(gmx_wallcycle_t wc)
249 if (wc->wcc != nullptr)
253 if (wc->wcc_all != nullptr)
257 if (wc->wcsc != nullptr)
264 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
267 wc->cycle_prev = cycle;
270 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
272 wc->wcc_all[wc->ewc_prev * ewcNR + ewc].n += 1;
273 wc->wcc_all[wc->ewc_prev * ewcNR + ewc].c += cycle - wc->cycle_prev;
278 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
280 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
282 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
284 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d", wc->count_depth);
286 wc->counterlist[wc->count_depth] = ewc;
290 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
294 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
296 if (wc->count_depth < 0)
298 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc],
301 if (wc->counterlist[wc->count_depth] != ewc)
303 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
304 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
309 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
321 MPI_Barrier(wc->mpi_comm_mygroup);
326 debug_start_check(wc, ewc);
329 cycle = gmx_cycles_read();
330 wc->wcc[ewc].start = cycle;
331 if (wc->wcc_all != nullptr)
336 wallcycle_all_start(wc, ewc, cycle);
338 else if (wc->wc_depth == 3)
340 wallcycle_all_stop(wc, ewc, cycle);
345 void wallcycle_increment_event_count(gmx_wallcycle_t wc, int ewc)
354 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
361 wallcycle_start(wc, ewc);
365 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
367 gmx_cycles_t cycle, last;
377 MPI_Barrier(wc->mpi_comm_mygroup);
382 debug_stop_check(wc, ewc);
385 /* When processes or threads migrate between cores, the cycle counting
386 * can get messed up if the cycle counter on different cores are not
387 * synchronized. When this happens we expect both large negative and
388 * positive cycle differences. We can detect negative cycle differences.
389 * Detecting too large positive counts if difficult, since count can be
390 * large, especially for ewcRUN. If we detect a negative count,
391 * we will not print the cycle accounting table.
393 cycle = gmx_cycles_read();
394 if (cycle >= wc->wcc[ewc].start)
396 last = cycle - wc->wcc[ewc].start;
401 wc->haveInvalidCount = TRUE;
403 wc->wcc[ewc].c += last;
410 wallcycle_all_stop(wc, ewc, cycle);
412 else if (wc->wc_depth == 2)
414 wallcycle_all_start(wc, ewc, cycle);
421 void wallcycle_get(gmx_wallcycle_t wc, int ewc, int* n, double* c)
424 *c = static_cast<double>(wc->wcc[ewc].c);
427 void wallcycle_reset_all(gmx_wallcycle_t wc)
436 for (i = 0; i < ewcNR; i++)
441 wc->haveInvalidCount = FALSE;
445 for (i = 0; i < ewcNR * ewcNR; i++)
447 wc->wcc_all[i].n = 0;
448 wc->wcc_all[i].c = 0;
453 for (i = 0; i < ewcsNR; i++)
461 static gmx_bool is_pme_counter(int ewc)
463 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
466 static gmx_bool is_pme_subcounter(int ewc)
468 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
471 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
472 static void subtract_cycles(wallcc_t* wcc, int ewc_main, int ewc_sub)
474 if (wcc[ewc_sub].n > 0)
476 if (wcc[ewc_main].c >= wcc[ewc_sub].c)
478 wcc[ewc_main].c -= wcc[ewc_sub].c;
482 /* Something is wrong with the cycle counting */
488 void wallcycle_scale_by_num_threads(gmx_wallcycle_t wc, bool isPmeRank, int nthreads_pp, int nthreads_pme)
495 for (int i = 0; i < ewcNR; i++)
497 if (is_pme_counter(i) || (i == ewcRUN && isPmeRank))
499 wc->wcc[i].c *= nthreads_pme;
503 for (int j = 0; j < ewcNR; j++)
505 wc->wcc_all[i * ewcNR + j].c *= nthreads_pme;
511 wc->wcc[i].c *= nthreads_pp;
515 for (int j = 0; j < ewcNR; j++)
517 wc->wcc_all[i * ewcNR + j].c *= nthreads_pp;
522 if (useCycleSubcounters && wc->wcsc && !isPmeRank)
524 for (int i = 0; i < ewcsNR; i++)
526 wc->wcsc[i].c *= nthreads_pp;
531 /* TODO Make an object for this function to return, containing some
532 * vectors of something like wallcc_t for the summed wcc, wcc_all and
533 * wcsc, AND the original wcc for rank 0.
535 * The GPU timing is reported only for rank 0, so we want to preserve
536 * the original wcycle on that rank. Rank 0 also reports the global
537 * counts before that, so needs something to contain the global data
538 * without over-writing the rank-0 data. The current implementation
539 * uses cycles_sum to manage this, which works OK now because wcsc and
540 * wcc_all are unused by the GPU reporting, but it is not satisfactory
541 * for the future. Also, there's no need for MPI_Allreduce, since
542 * only MASTERRANK uses any of the results. */
543 WallcycleCounts wallcycle_sum(const t_commrec* cr, gmx_wallcycle_t wc)
545 WallcycleCounts cycles_sum;
547 double cycles[ewcNR + ewcsNR];
549 double cycles_n[ewcNR + ewcsNR + 1];
556 /* Default construction of std::array of non-class T can leave
557 the values indeterminate, just like a C array, and icc
565 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
566 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
568 subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
570 if (cr->npmenodes == 0)
572 /* All nodes do PME (or no PME at all) */
573 subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
577 /* The are PME-only nodes */
578 if (wcc[ewcPMEMESH].n > 0)
580 /* This must be a PME only node, calculate the Wait + Comm. time */
581 GMX_ASSERT(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c,
582 "Total run ticks must be greater than PME-only ticks");
583 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
587 /* Store the cycles in a double buffer for summing */
588 for (i = 0; i < ewcNR; i++)
591 cycles_n[i] = static_cast<double>(wcc[i].n);
593 cycles[i] = static_cast<double>(wcc[i].c);
598 for (i = 0; i < ewcsNR; i++)
601 cycles_n[ewcNR + i] = static_cast<double>(wc->wcsc[i].n);
603 cycles[ewcNR + i] = static_cast<double>(wc->wcsc[i].c);
611 double buf[ewcNR + ewcsNR + 1];
613 // TODO this code is used only at the end of the run, so we
614 // can just do a simple reduce of haveInvalidCount in
615 // wallcycle_print, and avoid bugs
616 cycles_n[nsum] = (wc->haveInvalidCount ? 1 : 0);
617 // TODO Use MPI_Reduce
618 MPI_Allreduce(cycles_n, buf, nsum + 1, MPI_DOUBLE, MPI_MAX, cr->mpi_comm_mysim);
619 for (i = 0; i < ewcNR; i++)
621 wcc[i].n = gmx::roundToInt(buf[i]);
623 wc->haveInvalidCount = (buf[nsum] > 0);
626 for (i = 0; i < ewcsNR; i++)
628 wc->wcsc[i].n = gmx::roundToInt(buf[ewcNR + i]);
632 // TODO Use MPI_Reduce
633 MPI_Allreduce(cycles, cycles_sum.data(), nsum, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
635 if (wc->wcc_all != nullptr)
637 double *buf_all, *cyc_all;
639 snew(cyc_all, ewcNR * ewcNR);
640 snew(buf_all, ewcNR * ewcNR);
641 for (i = 0; i < ewcNR * ewcNR; i++)
643 cyc_all[i] = wc->wcc_all[i].c;
645 // TODO Use MPI_Reduce
646 MPI_Allreduce(cyc_all, buf_all, ewcNR * ewcNR, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
647 for (i = 0; i < ewcNR * ewcNR; i++)
649 wc->wcc_all[i].c = static_cast<gmx_cycles_t>(buf_all[i]);
658 for (i = 0; i < nsum; i++)
660 cycles_sum[i] = cycles[i];
668 print_cycles(FILE* fplog, double c2t, const char* name, int nnodes, int nthreads, int ncalls, double c_sum, double tot)
670 char nnodes_str[STRLEN];
671 char nthreads_str[STRLEN];
672 char ncalls_str[STRLEN];
674 double percentage = (tot > 0.) ? (100. * c_sum / tot) : 0.;
680 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
683 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
687 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
691 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
695 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
704 /* Convert the cycle count to wallclock time for this task */
707 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n", name, nnodes_str,
708 nthreads_str, ncalls_str, wallt, c_sum * 1e-9, percentage);
712 static void print_gputimes(FILE* fplog, const char* name, int n, double t, double tot_t)
719 snprintf(num, sizeof(num), "%10d", n);
720 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t / n);
725 sprintf(avg_perf, " ");
727 if (t != tot_t && tot_t > 0)
729 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n", name, num, t / 1000, avg_perf, 100 * t / tot_t);
733 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n", name, "", t / 1000, avg_perf, 100.0);
737 static void print_header(FILE* fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
739 int nrank_tot = nrank_pp + nrank_pme;
742 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
745 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
747 /* Don't report doing PP+PME, because we can't tell here if
748 * this is RF, etc. */
752 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
755 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
757 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
760 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
764 fprintf(fplog, "\n\n");
765 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
766 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
770 void wallcycle_print(FILE* fplog,
771 const gmx::MDLogger& mdlog,
778 const WallcycleCounts& cyc_sum,
779 const gmx_wallclock_gpu_nbnxn_t* gpu_nbnxn_t,
780 const gmx_wallclock_gpu_pme_t* gpu_pme_t)
782 double tot, tot_for_pp, tot_for_rest, tot_cpu_overlap, gpu_cpu_ratio;
783 double c2t, c2t_pp, c2t_pme = 0;
784 int i, j, npp, nth_tot;
787 "-----------------------------------------------------------------------------";
794 GMX_ASSERT(nth_pp > 0, "Number of particle-particle threads must be >0");
795 GMX_ASSERT(nth_pme > 0, "Number of PME threads must be >0");
796 GMX_ASSERT(nnodes > 0, "Number of nodes must be >0");
797 GMX_ASSERT(npme >= 0, "Number of PME nodes cannot be negative");
799 /* npme is the number of PME-only ranks used, and we always do PP work */
800 GMX_ASSERT(npp > 0, "Number of particle-particle nodes must be >0");
802 nth_tot = npp * nth_pp + npme * nth_pme;
804 /* When using PME-only nodes, the next line is valid for both
805 PP-only and PME-only nodes because they started ewcRUN at the
807 tot = cyc_sum[ewcRUN];
812 /* TODO This is heavy handed, but until someone reworks the
813 code so that it is provably robust with respect to
814 non-positive values for all possible timer and cycle
815 counters, there is less value gained from printing whatever
816 timing data might still be sensible for some non-Jenkins
817 run, than is lost from diagnosing Jenkins FP exceptions on
818 runs about whose execution time we don't care. */
819 GMX_LOG(mdlog.warning)
821 .appendTextFormatted(
822 "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a "
828 if (wc->haveInvalidCount)
830 GMX_LOG(mdlog.warning)
833 "NOTE: Detected invalid cycle counts, probably because threads moved "
834 "between CPU cores that do not have synchronized cycle counters. Will not "
835 "print the cycle accounting.");
840 /* Conversion factor from cycles to seconds */
841 c2t = realtime / tot;
842 c2t_pp = c2t * nth_tot / static_cast<double>(npp * nth_pp);
845 c2t_pme = c2t * nth_tot / static_cast<double>(npme * nth_pme);
852 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");
854 print_header(fplog, npp, nth_pp, npme, nth_pme);
856 fprintf(fplog, "%s\n", hline);
857 for (i = ewcPPDURINGPME + 1; i < ewcNR; i++)
859 if (is_pme_subcounter(i))
861 /* Do not count these at all */
863 else if (npme > 0 && is_pme_counter(i))
865 /* Print timing information for PME-only nodes, but add an
866 * asterisk so the reader of the table can know that the
867 * walltimes are not meant to add up. The asterisk still
868 * fits in the required maximum of 19 characters. */
870 snprintf(buffer, STRLEN, "%s *", wcn[i]);
871 print_cycles(fplog, c2t_pme, buffer, npme, nth_pme, wc->wcc[i].n, cyc_sum[i], tot);
875 /* Print timing information when it is for a PP or PP+PME
877 print_cycles(fplog, c2t_pp, wcn[i], npp, nth_pp, wc->wcc[i].n, cyc_sum[i], tot);
878 tot_for_pp += cyc_sum[i];
881 if (wc->wcc_all != nullptr)
883 for (i = 0; i < ewcNR; i++)
885 for (j = 0; j < ewcNR; j++)
887 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
888 print_cycles(fplog, c2t_pp, buf, npp, nth_pp, wc->wcc_all[i * ewcNR + j].n,
889 wc->wcc_all[i * ewcNR + j].c, tot);
893 tot_for_rest = tot * npp * nth_pp / static_cast<double>(nth_tot);
894 print_cycles(fplog, c2t_pp, "Rest", npp, nth_pp, -1, tot_for_rest - tot_for_pp, tot);
895 fprintf(fplog, "%s\n", hline);
896 print_cycles(fplog, c2t, "Total", npp, nth_pp, -1, tot, tot);
897 fprintf(fplog, "%s\n", hline);
902 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
903 " twice the total reported, but the cycle count total and %% are correct.\n"
908 if (wc->wcc[ewcPMEMESH].n > 0)
910 // A workaround to not print breakdown when no subcounters were recorded.
911 // TODO: figure out and record PME GPU counters (what to do with the waiting ones?)
912 std::vector<int> validPmeSubcounterIndices;
913 for (i = ewcPPDURINGPME + 1; i < ewcNR; i++)
915 if (is_pme_subcounter(i) && wc->wcc[i].n > 0)
917 validPmeSubcounterIndices.push_back(i);
921 if (!validPmeSubcounterIndices.empty())
923 fprintf(fplog, " Breakdown of PME mesh computation\n");
924 fprintf(fplog, "%s\n", hline);
925 for (auto i : validPmeSubcounterIndices)
927 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i], npme > 0 ? npme : npp,
928 nth_pme, wc->wcc[i].n, cyc_sum[i], tot);
930 fprintf(fplog, "%s\n", hline);
934 if (useCycleSubcounters && wc->wcsc)
936 fprintf(fplog, " Breakdown of PP computation\n");
937 fprintf(fplog, "%s\n", hline);
938 for (i = 0; i < ewcsNR; i++)
940 print_cycles(fplog, c2t_pp, wcsn[i], npp, nth_pp, wc->wcsc[i].n, cyc_sum[ewcNR + i], tot);
942 fprintf(fplog, "%s\n", hline);
945 /* print GPU timing summary */
946 double tot_gpu = 0.0;
949 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
951 tot_gpu += gpu_pme_t->timing[k].t;
956 const char* k_log_str[2][2] = { { "Nonbonded F kernel", "Nonbonded F+ene k." },
957 { "Nonbonded F+prune k.", "Nonbonded F+ene+prune k." } };
958 tot_gpu += gpu_nbnxn_t->pl_h2d_t + gpu_nbnxn_t->nb_h2d_t + gpu_nbnxn_t->nb_d2h_t;
960 /* add up the kernel timings */
961 for (i = 0; i < 2; i++)
963 for (j = 0; j < 2; j++)
965 tot_gpu += gpu_nbnxn_t->ktime[i][j].t;
968 tot_gpu += gpu_nbnxn_t->pruneTime.t;
970 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
971 if (wc->wcc[ewcPMEMESH].n > 0)
973 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
975 tot_cpu_overlap *= realtime * 1000 / tot; /* convert s to ms */
977 fprintf(fplog, "\n GPU timings\n%s\n", hline);
979 " Computing: Count Wall t (s) ms/step %c\n", '%');
980 fprintf(fplog, "%s\n", hline);
981 print_gputimes(fplog, "Pair list H2D", gpu_nbnxn_t->pl_h2d_c, gpu_nbnxn_t->pl_h2d_t, tot_gpu);
982 print_gputimes(fplog, "X / q H2D", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_h2d_t, tot_gpu);
984 for (i = 0; i < 2; i++)
986 for (j = 0; j < 2; j++)
988 if (gpu_nbnxn_t->ktime[i][j].c)
990 print_gputimes(fplog, k_log_str[i][j], gpu_nbnxn_t->ktime[i][j].c,
991 gpu_nbnxn_t->ktime[i][j].t, tot_gpu);
997 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
999 if (gpu_pme_t->timing[k].c)
1001 print_gputimes(fplog, PMEStageNames[k], gpu_pme_t->timing[k].c,
1002 gpu_pme_t->timing[k].t, tot_gpu);
1006 if (gpu_nbnxn_t->pruneTime.c)
1008 print_gputimes(fplog, "Pruning kernel", gpu_nbnxn_t->pruneTime.c,
1009 gpu_nbnxn_t->pruneTime.t, tot_gpu);
1011 print_gputimes(fplog, "F D2H", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_d2h_t, tot_gpu);
1012 fprintf(fplog, "%s\n", hline);
1013 print_gputimes(fplog, "Total ", gpu_nbnxn_t->nb_c, tot_gpu, tot_gpu);
1014 fprintf(fplog, "%s\n", hline);
1015 if (gpu_nbnxn_t->dynamicPruneTime.c)
1017 /* We print the dynamic pruning kernel timings after a separator
1018 * and avoid adding it to tot_gpu as this is not in the force
1019 * overlap. We print the fraction as relative to the rest.
1021 print_gputimes(fplog, "*Dynamic pruning", gpu_nbnxn_t->dynamicPruneTime.c,
1022 gpu_nbnxn_t->dynamicPruneTime.t, tot_gpu);
1023 fprintf(fplog, "%s\n", hline);
1025 gpu_cpu_ratio = tot_gpu / tot_cpu_overlap;
1026 if (gpu_nbnxn_t->nb_c > 0 && wc->wcc[ewcFORCE].n > 0)
1029 "\nAverage per-step force GPU/CPU evaluation time ratio: %.3f ms/%.3f ms = "
1031 tot_gpu / gpu_nbnxn_t->nb_c, tot_cpu_overlap / wc->wcc[ewcFORCE].n, gpu_cpu_ratio);
1034 /* only print notes related to CPU-GPU load balance with PME */
1035 if (wc->wcc[ewcPMEMESH].n > 0)
1037 fprintf(fplog, "For optimal resource utilization this ratio should be close to 1\n");
1039 /* print note if the imbalance is high with PME case in which
1040 * CPU-GPU load balancing is possible */
1041 if (gpu_cpu_ratio < 0.8 || gpu_cpu_ratio > 1.25)
1043 /* Only the sim master calls this function, so always print to stderr */
1044 if (gpu_cpu_ratio < 0.8)
1048 /* The user could have used -notunepme,
1049 * but we currently can't check that here.
1051 GMX_LOG(mdlog.warning)
1054 "NOTE: The CPU has >25% more load than the GPU. This "
1055 "imbalance wastes\n"
1056 " GPU resources. Maybe the domain decomposition "
1057 "limits the PME tuning.\n"
1058 " In that case, try setting the DD grid manually "
1059 "(-dd) or lowering -dds.");
1063 /* We should not end up here, unless the box is
1064 * too small for increasing the cut-off for PME tuning.
1066 GMX_LOG(mdlog.warning)
1069 "NOTE: The CPU has >25% more load than the GPU. This "
1070 "imbalance wastes\n"
1074 if (gpu_cpu_ratio > 1.25)
1076 GMX_LOG(mdlog.warning)
1079 "NOTE: The GPU has >25% more load than the CPU. This imbalance "
1089 GMX_LOG(mdlog.warning)
1092 "MPI_Barrier was called before each cycle start/stop\n"
1093 "call, so timings are not those of real runs.");
1096 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 && (cyc_sum[ewcDOMDEC] > tot * 0.1 || cyc_sum[ewcNS] > tot * 0.1))
1098 /* Only the sim master calls this function, so always print to stderr */
1099 if (wc->wcc[ewcDOMDEC].n == 0)
1101 GMX_LOG(mdlog.warning)
1103 .appendTextFormatted(
1104 "NOTE: %d %% of the run time was spent in pair search,\n"
1105 " you might want to increase nstlist (this has no effect on "
1107 gmx::roundToInt(100 * cyc_sum[ewcNS] / tot));
1111 GMX_LOG(mdlog.warning)
1113 .appendTextFormatted(
1114 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
1115 " %d %% of the run time was spent in pair search,\n"
1116 " you might want to increase nstlist (this has no effect on "
1118 gmx::roundToInt(100 * cyc_sum[ewcDOMDEC] / tot),
1119 gmx::roundToInt(100 * cyc_sum[ewcNS] / tot));
1123 if (cyc_sum[ewcMoveE] > tot * 0.05)
1125 GMX_LOG(mdlog.warning)
1127 .appendTextFormatted(
1128 "NOTE: %d %% of the run time was spent communicating energies,\n"
1129 " you might want to increase some nst* mdp options\n",
1130 gmx::roundToInt(100 * cyc_sum[ewcMoveE] / tot));
1134 extern int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
1141 return wc->reset_counters;
1144 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, int64_t reset_counters)
1151 wc->reset_counters = reset_counters;
1154 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
1156 if (useCycleSubcounters && wc != nullptr)
1158 wc->wcsc[ewcs].start = gmx_cycles_read();
1162 void wallcycle_sub_start_nocount(gmx_wallcycle_t wc, int ewcs)
1164 if (useCycleSubcounters && wc != nullptr)
1166 wallcycle_sub_start(wc, ewcs);
1171 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
1173 if (useCycleSubcounters && wc != nullptr)
1175 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;