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,2020, 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 "Wait GPU state copy",
136 "NB X/F buffer ops.",
150 static const char* wcsn[ewcsNR] = {
160 "NS search non-loc.",
164 "Listed buffer ops.",
166 "Nonbonded F kernel",
169 "Launch NB GPU tasks",
170 "Launch Bonded GPU tasks",
171 "Launch PME GPU tasks",
173 "Ewald F correction",
176 "Clear force buffer",
180 /* PME GPU timing events' names - correspond to the enum in the gpu_timing.h */
181 static const char* PMEStageNames[] = {
182 "PME spline", "PME spread", "PME spline + spread", "PME 3D-FFT r2c",
183 "PME solve", "PME 3D-FFT c2r", "PME gather",
186 gmx_bool wallcycle_have_counter()
188 return gmx_cycles_have_counter();
191 gmx_wallcycle_t wallcycle_init(FILE* fplog, int resetstep, t_commrec gmx_unused* cr)
196 if (!wallcycle_have_counter())
203 wc->haveInvalidCount = FALSE;
204 wc->wc_barrier = FALSE;
205 wc->wcc_all = nullptr;
208 wc->reset_counters = resetstep;
211 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != nullptr)
215 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
217 wc->wc_barrier = TRUE;
218 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
222 snew(wc->wcc, ewcNR);
223 if (getenv("GMX_CYCLE_ALL") != nullptr)
227 fprintf(fplog, "\nWill time all the code during the run\n\n");
229 snew(wc->wcc_all, ewcNR * ewcNR);
232 if (useCycleSubcounters)
234 snew(wc->wcsc, ewcsNR);
244 void wallcycle_destroy(gmx_wallcycle_t wc)
251 if (wc->wcc != nullptr)
255 if (wc->wcc_all != nullptr)
259 if (wc->wcsc != nullptr)
266 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
269 wc->cycle_prev = cycle;
272 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
274 wc->wcc_all[wc->ewc_prev * ewcNR + ewc].n += 1;
275 wc->wcc_all[wc->ewc_prev * ewcNR + ewc].c += cycle - wc->cycle_prev;
280 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
282 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
284 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
286 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d", wc->count_depth);
288 wc->counterlist[wc->count_depth] = ewc;
292 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
296 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
298 if (wc->count_depth < 0)
300 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc],
303 if (wc->counterlist[wc->count_depth] != ewc)
305 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
306 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
311 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
323 MPI_Barrier(wc->mpi_comm_mygroup);
328 debug_start_check(wc, ewc);
331 cycle = gmx_cycles_read();
332 wc->wcc[ewc].start = cycle;
333 if (wc->wcc_all != nullptr)
338 wallcycle_all_start(wc, ewc, cycle);
340 else if (wc->wc_depth == 3)
342 wallcycle_all_stop(wc, ewc, cycle);
347 void wallcycle_increment_event_count(gmx_wallcycle_t wc, int ewc)
356 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
363 wallcycle_start(wc, ewc);
367 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
369 gmx_cycles_t cycle, last;
379 MPI_Barrier(wc->mpi_comm_mygroup);
384 debug_stop_check(wc, ewc);
387 /* When processes or threads migrate between cores, the cycle counting
388 * can get messed up if the cycle counter on different cores are not
389 * synchronized. When this happens we expect both large negative and
390 * positive cycle differences. We can detect negative cycle differences.
391 * Detecting too large positive counts if difficult, since count can be
392 * large, especially for ewcRUN. If we detect a negative count,
393 * we will not print the cycle accounting table.
395 cycle = gmx_cycles_read();
396 if (cycle >= wc->wcc[ewc].start)
398 last = cycle - wc->wcc[ewc].start;
403 wc->haveInvalidCount = TRUE;
405 wc->wcc[ewc].c += last;
412 wallcycle_all_stop(wc, ewc, cycle);
414 else if (wc->wc_depth == 2)
416 wallcycle_all_start(wc, ewc, cycle);
423 void wallcycle_get(gmx_wallcycle_t wc, int ewc, int* n, double* c)
426 *c = static_cast<double>(wc->wcc[ewc].c);
429 void wallcycle_reset_all(gmx_wallcycle_t wc)
438 for (i = 0; i < ewcNR; i++)
443 wc->haveInvalidCount = FALSE;
447 for (i = 0; i < ewcNR * ewcNR; i++)
449 wc->wcc_all[i].n = 0;
450 wc->wcc_all[i].c = 0;
455 for (i = 0; i < ewcsNR; i++)
463 static gmx_bool is_pme_counter(int ewc)
465 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
468 static gmx_bool is_pme_subcounter(int ewc)
470 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
473 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
474 static void subtract_cycles(wallcc_t* wcc, int ewc_main, int ewc_sub)
476 if (wcc[ewc_sub].n > 0)
478 if (wcc[ewc_main].c >= wcc[ewc_sub].c)
480 wcc[ewc_main].c -= wcc[ewc_sub].c;
484 /* Something is wrong with the cycle counting */
490 void wallcycle_scale_by_num_threads(gmx_wallcycle_t wc, bool isPmeRank, int nthreads_pp, int nthreads_pme)
497 for (int i = 0; i < ewcNR; i++)
499 if (is_pme_counter(i) || (i == ewcRUN && isPmeRank))
501 wc->wcc[i].c *= nthreads_pme;
505 for (int j = 0; j < ewcNR; j++)
507 wc->wcc_all[i * ewcNR + j].c *= nthreads_pme;
513 wc->wcc[i].c *= nthreads_pp;
517 for (int j = 0; j < ewcNR; j++)
519 wc->wcc_all[i * ewcNR + j].c *= nthreads_pp;
524 if (useCycleSubcounters && wc->wcsc && !isPmeRank)
526 for (int i = 0; i < ewcsNR; i++)
528 wc->wcsc[i].c *= nthreads_pp;
533 /* TODO Make an object for this function to return, containing some
534 * vectors of something like wallcc_t for the summed wcc, wcc_all and
535 * wcsc, AND the original wcc for rank 0.
537 * The GPU timing is reported only for rank 0, so we want to preserve
538 * the original wcycle on that rank. Rank 0 also reports the global
539 * counts before that, so needs something to contain the global data
540 * without over-writing the rank-0 data. The current implementation
541 * uses cycles_sum to manage this, which works OK now because wcsc and
542 * wcc_all are unused by the GPU reporting, but it is not satisfactory
543 * for the future. Also, there's no need for MPI_Allreduce, since
544 * only MASTERRANK uses any of the results. */
545 WallcycleCounts wallcycle_sum(const t_commrec* cr, gmx_wallcycle_t wc)
547 WallcycleCounts cycles_sum;
549 double cycles[ewcNR + ewcsNR];
551 double cycles_n[ewcNR + ewcsNR + 1];
558 /* Default construction of std::array of non-class T can leave
559 the values indeterminate, just like a C array, and icc
567 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
568 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
570 subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
572 if (cr->npmenodes == 0)
574 /* All nodes do PME (or no PME at all) */
575 subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
579 /* The are PME-only nodes */
580 if (wcc[ewcPMEMESH].n > 0)
582 /* This must be a PME only node, calculate the Wait + Comm. time */
583 GMX_ASSERT(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c,
584 "Total run ticks must be greater than PME-only ticks");
585 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
589 /* Store the cycles in a double buffer for summing */
590 for (i = 0; i < ewcNR; i++)
593 cycles_n[i] = static_cast<double>(wcc[i].n);
595 cycles[i] = static_cast<double>(wcc[i].c);
600 for (i = 0; i < ewcsNR; i++)
603 cycles_n[ewcNR + i] = static_cast<double>(wc->wcsc[i].n);
605 cycles[ewcNR + i] = static_cast<double>(wc->wcsc[i].c);
613 double buf[ewcNR + ewcsNR + 1];
615 // TODO this code is used only at the end of the run, so we
616 // can just do a simple reduce of haveInvalidCount in
617 // wallcycle_print, and avoid bugs
618 cycles_n[nsum] = (wc->haveInvalidCount ? 1 : 0);
619 // TODO Use MPI_Reduce
620 MPI_Allreduce(cycles_n, buf, nsum + 1, MPI_DOUBLE, MPI_MAX, cr->mpi_comm_mysim);
621 for (i = 0; i < ewcNR; i++)
623 wcc[i].n = gmx::roundToInt(buf[i]);
625 wc->haveInvalidCount = (buf[nsum] > 0);
628 for (i = 0; i < ewcsNR; i++)
630 wc->wcsc[i].n = gmx::roundToInt(buf[ewcNR + i]);
634 // TODO Use MPI_Reduce
635 MPI_Allreduce(cycles, cycles_sum.data(), nsum, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
637 if (wc->wcc_all != nullptr)
639 double *buf_all, *cyc_all;
641 snew(cyc_all, ewcNR * ewcNR);
642 snew(buf_all, ewcNR * ewcNR);
643 for (i = 0; i < ewcNR * ewcNR; i++)
645 cyc_all[i] = wc->wcc_all[i].c;
647 // TODO Use MPI_Reduce
648 MPI_Allreduce(cyc_all, buf_all, ewcNR * ewcNR, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
649 for (i = 0; i < ewcNR * ewcNR; i++)
651 wc->wcc_all[i].c = static_cast<gmx_cycles_t>(buf_all[i]);
660 for (i = 0; i < nsum; i++)
662 cycles_sum[i] = cycles[i];
670 print_cycles(FILE* fplog, double c2t, const char* name, int nnodes, int nthreads, int ncalls, double c_sum, double tot)
672 char nnodes_str[STRLEN];
673 char nthreads_str[STRLEN];
674 char ncalls_str[STRLEN];
676 double percentage = (tot > 0.) ? (100. * c_sum / tot) : 0.;
682 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
685 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
689 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
693 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
697 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
706 /* Convert the cycle count to wallclock time for this task */
709 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n", name, nnodes_str,
710 nthreads_str, ncalls_str, wallt, c_sum * 1e-9, percentage);
714 static void print_gputimes(FILE* fplog, const char* name, int n, double t, double tot_t)
721 snprintf(num, sizeof(num), "%10d", n);
722 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t / n);
727 sprintf(avg_perf, " ");
729 if (t != tot_t && tot_t > 0)
731 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n", name, num, t / 1000, avg_perf, 100 * t / tot_t);
735 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n", name, "", t / 1000, avg_perf, 100.0);
739 static void print_header(FILE* fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
741 int nrank_tot = nrank_pp + nrank_pme;
744 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
747 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
749 /* Don't report doing PP+PME, because we can't tell here if
750 * this is RF, etc. */
754 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
757 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
759 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
762 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
766 fprintf(fplog, "\n\n");
767 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
768 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
772 void wallcycle_print(FILE* fplog,
773 const gmx::MDLogger& mdlog,
780 const WallcycleCounts& cyc_sum,
781 const gmx_wallclock_gpu_nbnxn_t* gpu_nbnxn_t,
782 const gmx_wallclock_gpu_pme_t* gpu_pme_t)
784 double tot, tot_for_pp, tot_for_rest, tot_cpu_overlap, gpu_cpu_ratio;
785 double c2t, c2t_pp, c2t_pme = 0;
786 int i, j, npp, nth_tot;
789 "-----------------------------------------------------------------------------";
796 GMX_ASSERT(nth_pp > 0, "Number of particle-particle threads must be >0");
797 GMX_ASSERT(nth_pme > 0, "Number of PME threads must be >0");
798 GMX_ASSERT(nnodes > 0, "Number of nodes must be >0");
799 GMX_ASSERT(npme >= 0, "Number of PME nodes cannot be negative");
801 /* npme is the number of PME-only ranks used, and we always do PP work */
802 GMX_ASSERT(npp > 0, "Number of particle-particle nodes must be >0");
804 nth_tot = npp * nth_pp + npme * nth_pme;
806 /* When using PME-only nodes, the next line is valid for both
807 PP-only and PME-only nodes because they started ewcRUN at the
809 tot = cyc_sum[ewcRUN];
814 /* TODO This is heavy handed, but until someone reworks the
815 code so that it is provably robust with respect to
816 non-positive values for all possible timer and cycle
817 counters, there is less value gained from printing whatever
818 timing data might still be sensible for some non-Jenkins
819 run, than is lost from diagnosing Jenkins FP exceptions on
820 runs about whose execution time we don't care. */
821 GMX_LOG(mdlog.warning)
823 .appendTextFormatted(
824 "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a "
830 if (wc->haveInvalidCount)
832 GMX_LOG(mdlog.warning)
835 "NOTE: Detected invalid cycle counts, probably because threads moved "
836 "between CPU cores that do not have synchronized cycle counters. Will not "
837 "print the cycle accounting.");
842 /* Conversion factor from cycles to seconds */
843 c2t = realtime / tot;
844 c2t_pp = c2t * nth_tot / static_cast<double>(npp * nth_pp);
847 c2t_pme = c2t * nth_tot / static_cast<double>(npme * nth_pme);
854 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");
856 print_header(fplog, npp, nth_pp, npme, nth_pme);
858 fprintf(fplog, "%s\n", hline);
859 for (i = ewcPPDURINGPME + 1; i < ewcNR; i++)
861 if (is_pme_subcounter(i))
863 /* Do not count these at all */
865 else if (npme > 0 && is_pme_counter(i))
867 /* Print timing information for PME-only nodes, but add an
868 * asterisk so the reader of the table can know that the
869 * walltimes are not meant to add up. The asterisk still
870 * fits in the required maximum of 19 characters. */
872 snprintf(buffer, STRLEN, "%s *", wcn[i]);
873 print_cycles(fplog, c2t_pme, buffer, npme, nth_pme, wc->wcc[i].n, cyc_sum[i], tot);
877 /* Print timing information when it is for a PP or PP+PME
879 print_cycles(fplog, c2t_pp, wcn[i], npp, nth_pp, wc->wcc[i].n, cyc_sum[i], tot);
880 tot_for_pp += cyc_sum[i];
883 if (wc->wcc_all != nullptr)
885 for (i = 0; i < ewcNR; i++)
887 for (j = 0; j < ewcNR; j++)
889 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
890 print_cycles(fplog, c2t_pp, buf, npp, nth_pp, wc->wcc_all[i * ewcNR + j].n,
891 wc->wcc_all[i * ewcNR + j].c, tot);
895 tot_for_rest = tot * npp * nth_pp / static_cast<double>(nth_tot);
896 print_cycles(fplog, c2t_pp, "Rest", npp, nth_pp, -1, tot_for_rest - tot_for_pp, tot);
897 fprintf(fplog, "%s\n", hline);
898 print_cycles(fplog, c2t, "Total", npp, nth_pp, -1, tot, tot);
899 fprintf(fplog, "%s\n", hline);
904 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
905 " twice the total reported, but the cycle count total and %% are correct.\n"
910 if (wc->wcc[ewcPMEMESH].n > 0)
912 // A workaround to not print breakdown when no subcounters were recorded.
913 // TODO: figure out and record PME GPU counters (what to do with the waiting ones?)
914 std::vector<int> validPmeSubcounterIndices;
915 for (i = ewcPPDURINGPME + 1; i < ewcNR; i++)
917 if (is_pme_subcounter(i) && wc->wcc[i].n > 0)
919 validPmeSubcounterIndices.push_back(i);
923 if (!validPmeSubcounterIndices.empty())
925 fprintf(fplog, " Breakdown of PME mesh computation\n");
926 fprintf(fplog, "%s\n", hline);
927 for (auto i : validPmeSubcounterIndices)
929 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i], npme > 0 ? npme : npp,
930 nth_pme, wc->wcc[i].n, cyc_sum[i], tot);
932 fprintf(fplog, "%s\n", hline);
936 if (useCycleSubcounters && wc->wcsc)
938 fprintf(fplog, " Breakdown of PP computation\n");
939 fprintf(fplog, "%s\n", hline);
940 for (i = 0; i < ewcsNR; i++)
942 print_cycles(fplog, c2t_pp, wcsn[i], npp, nth_pp, wc->wcsc[i].n, cyc_sum[ewcNR + i], tot);
944 fprintf(fplog, "%s\n", hline);
947 /* print GPU timing summary */
948 double tot_gpu = 0.0;
951 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
953 tot_gpu += gpu_pme_t->timing[k].t;
958 const char* k_log_str[2][2] = { { "Nonbonded F kernel", "Nonbonded F+ene k." },
959 { "Nonbonded F+prune k.", "Nonbonded F+ene+prune k." } };
960 tot_gpu += gpu_nbnxn_t->pl_h2d_t + gpu_nbnxn_t->nb_h2d_t + gpu_nbnxn_t->nb_d2h_t;
962 /* add up the kernel timings */
963 for (i = 0; i < 2; i++)
965 for (j = 0; j < 2; j++)
967 tot_gpu += gpu_nbnxn_t->ktime[i][j].t;
970 tot_gpu += gpu_nbnxn_t->pruneTime.t;
972 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
973 if (wc->wcc[ewcPMEMESH].n > 0)
975 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
977 tot_cpu_overlap *= realtime * 1000 / tot; /* convert s to ms */
979 fprintf(fplog, "\n GPU timings\n%s\n", hline);
981 " Computing: Count Wall t (s) ms/step %c\n", '%');
982 fprintf(fplog, "%s\n", hline);
983 print_gputimes(fplog, "Pair list H2D", gpu_nbnxn_t->pl_h2d_c, gpu_nbnxn_t->pl_h2d_t, tot_gpu);
984 print_gputimes(fplog, "X / q H2D", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_h2d_t, tot_gpu);
986 for (i = 0; i < 2; i++)
988 for (j = 0; j < 2; j++)
990 if (gpu_nbnxn_t->ktime[i][j].c)
992 print_gputimes(fplog, k_log_str[i][j], gpu_nbnxn_t->ktime[i][j].c,
993 gpu_nbnxn_t->ktime[i][j].t, tot_gpu);
999 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
1001 if (gpu_pme_t->timing[k].c)
1003 print_gputimes(fplog, PMEStageNames[k], gpu_pme_t->timing[k].c,
1004 gpu_pme_t->timing[k].t, tot_gpu);
1008 if (gpu_nbnxn_t->pruneTime.c)
1010 print_gputimes(fplog, "Pruning kernel", gpu_nbnxn_t->pruneTime.c,
1011 gpu_nbnxn_t->pruneTime.t, tot_gpu);
1013 print_gputimes(fplog, "F D2H", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_d2h_t, tot_gpu);
1014 fprintf(fplog, "%s\n", hline);
1015 print_gputimes(fplog, "Total ", gpu_nbnxn_t->nb_c, tot_gpu, tot_gpu);
1016 fprintf(fplog, "%s\n", hline);
1017 if (gpu_nbnxn_t->dynamicPruneTime.c)
1019 /* We print the dynamic pruning kernel timings after a separator
1020 * and avoid adding it to tot_gpu as this is not in the force
1021 * overlap. We print the fraction as relative to the rest.
1023 print_gputimes(fplog, "*Dynamic pruning", gpu_nbnxn_t->dynamicPruneTime.c,
1024 gpu_nbnxn_t->dynamicPruneTime.t, tot_gpu);
1025 fprintf(fplog, "%s\n", hline);
1027 gpu_cpu_ratio = tot_gpu / tot_cpu_overlap;
1028 if (gpu_nbnxn_t->nb_c > 0 && wc->wcc[ewcFORCE].n > 0)
1031 "\nAverage per-step force GPU/CPU evaluation time ratio: %.3f ms/%.3f ms = "
1033 tot_gpu / gpu_nbnxn_t->nb_c, tot_cpu_overlap / wc->wcc[ewcFORCE].n, gpu_cpu_ratio);
1036 /* only print notes related to CPU-GPU load balance with PME */
1037 if (wc->wcc[ewcPMEMESH].n > 0)
1039 fprintf(fplog, "For optimal resource utilization this ratio should be close to 1\n");
1041 /* print note if the imbalance is high with PME case in which
1042 * CPU-GPU load balancing is possible */
1043 if (gpu_cpu_ratio < 0.8 || gpu_cpu_ratio > 1.25)
1045 /* Only the sim master calls this function, so always print to stderr */
1046 if (gpu_cpu_ratio < 0.8)
1050 /* The user could have used -notunepme,
1051 * but we currently can't check that here.
1053 GMX_LOG(mdlog.warning)
1056 "NOTE: The CPU has >25% more load than the GPU. This "
1057 "imbalance wastes\n"
1058 " GPU resources. Maybe the domain decomposition "
1059 "limits the PME tuning.\n"
1060 " In that case, try setting the DD grid manually "
1061 "(-dd) or lowering -dds.");
1065 /* We should not end up here, unless the box is
1066 * too small for increasing the cut-off for PME tuning.
1068 GMX_LOG(mdlog.warning)
1071 "NOTE: The CPU has >25% more load than the GPU. This "
1072 "imbalance wastes\n"
1076 if (gpu_cpu_ratio > 1.25)
1078 GMX_LOG(mdlog.warning)
1081 "NOTE: The GPU has >25% more load than the CPU. This imbalance "
1091 GMX_LOG(mdlog.warning)
1094 "MPI_Barrier was called before each cycle start/stop\n"
1095 "call, so timings are not those of real runs.");
1098 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 && (cyc_sum[ewcDOMDEC] > tot * 0.1 || cyc_sum[ewcNS] > tot * 0.1))
1100 /* Only the sim master calls this function, so always print to stderr */
1101 if (wc->wcc[ewcDOMDEC].n == 0)
1103 GMX_LOG(mdlog.warning)
1105 .appendTextFormatted(
1106 "NOTE: %d %% of the run time was spent in pair search,\n"
1107 " you might want to increase nstlist (this has no effect on "
1109 gmx::roundToInt(100 * cyc_sum[ewcNS] / tot));
1113 GMX_LOG(mdlog.warning)
1115 .appendTextFormatted(
1116 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
1117 " %d %% of the run time was spent in pair search,\n"
1118 " you might want to increase nstlist (this has no effect on "
1120 gmx::roundToInt(100 * cyc_sum[ewcDOMDEC] / tot),
1121 gmx::roundToInt(100 * cyc_sum[ewcNS] / tot));
1125 if (cyc_sum[ewcMoveE] > tot * 0.05)
1127 GMX_LOG(mdlog.warning)
1129 .appendTextFormatted(
1130 "NOTE: %d %% of the run time was spent communicating energies,\n"
1131 " you might want to increase some nst* mdp options\n",
1132 gmx::roundToInt(100 * cyc_sum[ewcMoveE] / tot));
1136 extern int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
1143 return wc->reset_counters;
1146 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, int64_t reset_counters)
1153 wc->reset_counters = reset_counters;
1156 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
1158 if (useCycleSubcounters && wc != nullptr)
1160 wc->wcsc[ewcs].start = gmx_cycles_read();
1164 void wallcycle_sub_start_nocount(gmx_wallcycle_t wc, int ewcs)
1166 if (useCycleSubcounters && wc != nullptr)
1168 wallcycle_sub_start(wc, ewcs);
1173 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
1175 if (useCycleSubcounters && wc != nullptr)
1177 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;