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] =
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.",
127 "Listed buffer ops.",
130 "Launch NB GPU tasks",
131 "Launch Bonded GPU tasks",
132 "Launch PME GPU tasks",
133 "Ewald F correction",
136 "Clear force buffer",
140 /* PME GPU timing events' names - correspond to the enum in the gpu_timing.h */
141 static const char *PMEStageNames[] =
145 "PME spline + spread",
152 gmx_bool wallcycle_have_counter()
154 return gmx_cycles_have_counter();
157 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec gmx_unused *cr)
162 if (!wallcycle_have_counter())
169 wc->haveInvalidCount = FALSE;
170 wc->wc_barrier = FALSE;
171 wc->wcc_all = nullptr;
174 wc->reset_counters = resetstep;
177 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != nullptr)
181 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
183 wc->wc_barrier = TRUE;
184 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
188 snew(wc->wcc, ewcNR);
189 if (getenv("GMX_CYCLE_ALL") != nullptr)
193 fprintf(fplog, "\nWill time all the code during the run\n\n");
195 snew(wc->wcc_all, ewcNR*ewcNR);
198 if (useCycleSubcounters)
200 snew(wc->wcsc, ewcsNR);
210 void wallcycle_destroy(gmx_wallcycle_t wc)
217 if (wc->wcc != nullptr)
221 if (wc->wcc_all != nullptr)
225 if (wc->wcsc != nullptr)
232 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
235 wc->cycle_prev = cycle;
238 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
240 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
241 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
246 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
248 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
250 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
252 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
255 wc->counterlist[wc->count_depth] = ewc;
259 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
263 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
265 if (wc->count_depth < 0)
267 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
269 if (wc->counterlist[wc->count_depth] != ewc)
271 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
272 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
277 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
289 MPI_Barrier(wc->mpi_comm_mygroup);
294 debug_start_check(wc, ewc);
297 cycle = gmx_cycles_read();
298 wc->wcc[ewc].start = cycle;
299 if (wc->wcc_all != nullptr)
304 wallcycle_all_start(wc, ewc, cycle);
306 else if (wc->wc_depth == 3)
308 wallcycle_all_stop(wc, ewc, cycle);
313 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
320 wallcycle_start(wc, ewc);
324 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
326 gmx_cycles_t cycle, last;
336 MPI_Barrier(wc->mpi_comm_mygroup);
341 debug_stop_check(wc, ewc);
344 /* When processes or threads migrate between cores, the cycle counting
345 * can get messed up if the cycle counter on different cores are not
346 * synchronized. When this happens we expect both large negative and
347 * positive cycle differences. We can detect negative cycle differences.
348 * Detecting too large positive counts if difficult, since count can be
349 * large, especially for ewcRUN. If we detect a negative count,
350 * we will not print the cycle accounting table.
352 cycle = gmx_cycles_read();
353 if (cycle >= wc->wcc[ewc].start)
355 last = cycle - wc->wcc[ewc].start;
360 wc->haveInvalidCount = TRUE;
362 wc->wcc[ewc].c += last;
369 wallcycle_all_stop(wc, ewc, cycle);
371 else if (wc->wc_depth == 2)
373 wallcycle_all_start(wc, ewc, cycle);
380 void wallcycle_get(gmx_wallcycle_t wc, int ewc, int *n, double *c)
383 *c = static_cast<double>(wc->wcc[ewc].c);
386 void wallcycle_reset_all(gmx_wallcycle_t wc)
395 for (i = 0; i < ewcNR; i++)
400 wc->haveInvalidCount = FALSE;
404 for (i = 0; i < ewcNR*ewcNR; i++)
406 wc->wcc_all[i].n = 0;
407 wc->wcc_all[i].c = 0;
412 for (i = 0; i < ewcsNR; i++)
420 static gmx_bool is_pme_counter(int ewc)
422 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
425 static gmx_bool is_pme_subcounter(int ewc)
427 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
430 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
431 static void subtract_cycles(wallcc_t *wcc, int ewc_main, int ewc_sub)
433 if (wcc[ewc_sub].n > 0)
435 if (wcc[ewc_main].c >= wcc[ewc_sub].c)
437 wcc[ewc_main].c -= wcc[ewc_sub].c;
441 /* Something is wrong with the cycle counting */
447 void wallcycle_scale_by_num_threads(gmx_wallcycle_t wc, bool isPmeRank, int nthreads_pp, int nthreads_pme)
454 for (int i = 0; i < ewcNR; i++)
456 if (is_pme_counter(i) || (i == ewcRUN && isPmeRank))
458 wc->wcc[i].c *= nthreads_pme;
462 for (int j = 0; j < ewcNR; j++)
464 wc->wcc_all[i*ewcNR+j].c *= nthreads_pme;
470 wc->wcc[i].c *= nthreads_pp;
474 for (int j = 0; j < ewcNR; j++)
476 wc->wcc_all[i*ewcNR+j].c *= nthreads_pp;
481 if (useCycleSubcounters && wc->wcsc && !isPmeRank)
483 for (int i = 0; i < ewcsNR; i++)
485 wc->wcsc[i].c *= nthreads_pp;
490 /* TODO Make an object for this function to return, containing some
491 * vectors of something like wallcc_t for the summed wcc, wcc_all and
492 * wcsc, AND the original wcc for rank 0.
494 * The GPU timing is reported only for rank 0, so we want to preserve
495 * the original wcycle on that rank. Rank 0 also reports the global
496 * counts before that, so needs something to contain the global data
497 * without over-writing the rank-0 data. The current implementation
498 * uses cycles_sum to manage this, which works OK now because wcsc and
499 * wcc_all are unused by the GPU reporting, but it is not satisfactory
500 * for the future. Also, there's no need for MPI_Allreduce, since
501 * only MASTERRANK uses any of the results. */
502 WallcycleCounts wallcycle_sum(const t_commrec *cr, gmx_wallcycle_t wc)
504 WallcycleCounts cycles_sum;
506 double cycles[ewcNR+ewcsNR];
508 double cycles_n[ewcNR+ewcsNR+1];
515 /* Default construction of std::array of non-class T can leave
516 the values indeterminate, just like a C array, and icc
524 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
525 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
527 subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
529 if (cr->npmenodes == 0)
531 /* All nodes do PME (or no PME at all) */
532 subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
536 /* The are PME-only nodes */
537 if (wcc[ewcPMEMESH].n > 0)
539 /* This must be a PME only node, calculate the Wait + Comm. time */
540 GMX_ASSERT(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c, "Total run ticks must be greater than PME-only ticks");
541 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
545 /* Store the cycles in a double buffer for summing */
546 for (i = 0; i < ewcNR; i++)
549 cycles_n[i] = static_cast<double>(wcc[i].n);
551 cycles[i] = static_cast<double>(wcc[i].c);
556 for (i = 0; i < ewcsNR; i++)
559 cycles_n[ewcNR+i] = static_cast<double>(wc->wcsc[i].n);
561 cycles[ewcNR+i] = static_cast<double>(wc->wcsc[i].c);
569 double buf[ewcNR+ewcsNR+1];
571 // TODO this code is used only at the end of the run, so we
572 // can just do a simple reduce of haveInvalidCount in
573 // wallcycle_print, and avoid bugs
574 cycles_n[nsum] = (wc->haveInvalidCount ? 1 : 0);
575 // TODO Use MPI_Reduce
576 MPI_Allreduce(cycles_n, buf, nsum + 1, MPI_DOUBLE, MPI_MAX,
578 for (i = 0; i < ewcNR; i++)
580 wcc[i].n = gmx::roundToInt(buf[i]);
582 wc->haveInvalidCount = (buf[nsum] > 0);
585 for (i = 0; i < ewcsNR; i++)
587 wc->wcsc[i].n = gmx::roundToInt(buf[ewcNR+i]);
591 // TODO Use MPI_Reduce
592 MPI_Allreduce(cycles, cycles_sum.data(), nsum, MPI_DOUBLE, MPI_SUM,
595 if (wc->wcc_all != nullptr)
597 double *buf_all, *cyc_all;
599 snew(cyc_all, ewcNR*ewcNR);
600 snew(buf_all, ewcNR*ewcNR);
601 for (i = 0; i < ewcNR*ewcNR; i++)
603 cyc_all[i] = wc->wcc_all[i].c;
605 // TODO Use MPI_Reduce
606 MPI_Allreduce(cyc_all, buf_all, ewcNR*ewcNR, MPI_DOUBLE, MPI_SUM,
608 for (i = 0; i < ewcNR*ewcNR; i++)
610 wc->wcc_all[i].c = static_cast<gmx_cycles_t>(buf_all[i]);
619 for (i = 0; i < nsum; i++)
621 cycles_sum[i] = cycles[i];
628 static void print_cycles(FILE *fplog, double c2t, const char *name,
629 int nnodes, int nthreads,
630 int ncalls, double c_sum, double tot)
632 char nnodes_str[STRLEN];
633 char nthreads_str[STRLEN];
634 char ncalls_str[STRLEN];
636 double percentage = (tot > 0.) ? (100. * c_sum / tot) : 0.;
642 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
645 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
649 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
653 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
657 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
666 /* Convert the cycle count to wallclock time for this task */
669 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
670 name, nnodes_str, nthreads_str, ncalls_str, wallt,
671 c_sum*1e-9, percentage);
675 static void print_gputimes(FILE *fplog, const char *name,
676 int n, double t, double tot_t)
683 snprintf(num, sizeof(num), "%10d", n);
684 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
689 sprintf(avg_perf, " ");
691 if (t != tot_t && tot_t > 0)
693 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
694 name, num, t/1000, avg_perf, 100 * t/tot_t);
698 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
699 name, "", t/1000, avg_perf, 100.0);
703 static void print_header(FILE *fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
705 int nrank_tot = nrank_pp + nrank_pme;
708 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
711 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
713 /* Don't report doing PP+PME, because we can't tell here if
714 * this is RF, etc. */
718 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
721 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
723 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
726 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
730 fprintf(fplog, "\n\n");
731 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
732 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
736 void wallcycle_print(FILE *fplog, const gmx::MDLogger &mdlog, int nnodes, int npme,
737 int nth_pp, int nth_pme, double realtime,
738 gmx_wallcycle_t wc, const WallcycleCounts &cyc_sum,
739 const gmx_wallclock_gpu_nbnxn_t *gpu_nbnxn_t,
740 const gmx_wallclock_gpu_pme_t *gpu_pme_t)
742 double tot, tot_for_pp, tot_for_rest, tot_cpu_overlap, gpu_cpu_ratio;
743 double c2t, c2t_pp, c2t_pme = 0;
744 int i, j, npp, nth_tot;
746 const char *hline = "-----------------------------------------------------------------------------";
753 GMX_ASSERT(nth_pp > 0, "Number of particle-particle threads must be >0");
754 GMX_ASSERT(nth_pme > 0, "Number of PME threads must be >0");
755 GMX_ASSERT(nnodes > 0, "Number of nodes must be >0");
756 GMX_ASSERT(npme >= 0, "Number of PME nodes cannot be negative");
758 /* npme is the number of PME-only ranks used, and we always do PP work */
759 GMX_ASSERT(npp > 0, "Number of particle-particle nodes must be >0");
761 nth_tot = npp*nth_pp + npme*nth_pme;
763 /* When using PME-only nodes, the next line is valid for both
764 PP-only and PME-only nodes because they started ewcRUN at the
766 tot = cyc_sum[ewcRUN];
771 /* TODO This is heavy handed, but until someone reworks the
772 code so that it is provably robust with respect to
773 non-positive values for all possible timer and cycle
774 counters, there is less value gained from printing whatever
775 timing data might still be sensible for some non-Jenkins
776 run, than is lost from diagnosing Jenkins FP exceptions on
777 runs about whose execution time we don't care. */
778 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
779 "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a time accounting",
784 if (wc->haveInvalidCount)
786 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.");
791 /* Conversion factor from cycles to seconds */
793 c2t_pp = c2t * nth_tot / static_cast<double>(npp*nth_pp);
796 c2t_pme = c2t * nth_tot / static_cast<double>(npme*nth_pme);
803 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");
805 print_header(fplog, npp, nth_pp, npme, nth_pme);
807 fprintf(fplog, "%s\n", hline);
808 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
810 if (is_pme_subcounter(i))
812 /* Do not count these at all */
814 else if (npme > 0 && is_pme_counter(i))
816 /* Print timing information for PME-only nodes, but add an
817 * asterisk so the reader of the table can know that the
818 * walltimes are not meant to add up. The asterisk still
819 * fits in the required maximum of 19 characters. */
821 snprintf(buffer, STRLEN, "%s *", wcn[i]);
822 print_cycles(fplog, c2t_pme, buffer,
824 wc->wcc[i].n, cyc_sum[i], tot);
828 /* Print timing information when it is for a PP or PP+PME
830 print_cycles(fplog, c2t_pp, wcn[i],
832 wc->wcc[i].n, cyc_sum[i], tot);
833 tot_for_pp += cyc_sum[i];
836 if (wc->wcc_all != nullptr)
838 for (i = 0; i < ewcNR; i++)
840 for (j = 0; j < ewcNR; j++)
842 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
843 print_cycles(fplog, c2t_pp, buf,
845 wc->wcc_all[i*ewcNR+j].n,
846 wc->wcc_all[i*ewcNR+j].c,
851 tot_for_rest = tot * npp * nth_pp / static_cast<double>(nth_tot);
852 print_cycles(fplog, c2t_pp, "Rest",
854 -1, tot_for_rest - tot_for_pp, tot);
855 fprintf(fplog, "%s\n", hline);
856 print_cycles(fplog, c2t, "Total",
859 fprintf(fplog, "%s\n", hline);
864 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
865 " twice the total reported, but the cycle count total and %% are correct.\n"
869 if (wc->wcc[ewcPMEMESH].n > 0)
871 // A workaround to not print breakdown when no subcounters were recorded.
872 // TODO: figure out and record PME GPU counters (what to do with the waiting ones?)
873 std::vector<int> validPmeSubcounterIndices;
874 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
876 if (is_pme_subcounter(i) && wc->wcc[i].n > 0)
878 validPmeSubcounterIndices.push_back(i);
882 if (!validPmeSubcounterIndices.empty())
884 fprintf(fplog, " Breakdown of PME mesh computation\n");
885 fprintf(fplog, "%s\n", hline);
886 for (auto i : validPmeSubcounterIndices)
888 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i],
889 npme > 0 ? npme : npp, nth_pme,
890 wc->wcc[i].n, cyc_sum[i], tot);
892 fprintf(fplog, "%s\n", hline);
896 if (useCycleSubcounters && wc->wcsc)
898 fprintf(fplog, " Breakdown of PP computation\n");
899 fprintf(fplog, "%s\n", hline);
900 for (i = 0; i < ewcsNR; i++)
902 print_cycles(fplog, c2t_pp, wcsn[i],
904 wc->wcsc[i].n, cyc_sum[ewcNR+i], tot);
906 fprintf(fplog, "%s\n", hline);
909 /* print GPU timing summary */
910 double tot_gpu = 0.0;
913 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
915 tot_gpu += gpu_pme_t->timing[k].t;
920 const char *k_log_str[2][2] = {
921 {"Nonbonded F kernel", "Nonbonded F+ene k."},
922 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
924 tot_gpu += gpu_nbnxn_t->pl_h2d_t + gpu_nbnxn_t->nb_h2d_t + gpu_nbnxn_t->nb_d2h_t;
926 /* add up the kernel timings */
927 for (i = 0; i < 2; i++)
929 for (j = 0; j < 2; j++)
931 tot_gpu += gpu_nbnxn_t->ktime[i][j].t;
934 tot_gpu += gpu_nbnxn_t->pruneTime.t;
936 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
937 if (wc->wcc[ewcPMEMESH].n > 0)
939 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
941 tot_cpu_overlap *= realtime*1000/tot; /* convert s to ms */
943 fprintf(fplog, "\n GPU timings\n%s\n", hline);
944 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
945 fprintf(fplog, "%s\n", hline);
946 print_gputimes(fplog, "Pair list H2D",
947 gpu_nbnxn_t->pl_h2d_c, gpu_nbnxn_t->pl_h2d_t, tot_gpu);
948 print_gputimes(fplog, "X / q H2D",
949 gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_h2d_t, tot_gpu);
951 for (i = 0; i < 2; i++)
953 for (j = 0; j < 2; j++)
955 if (gpu_nbnxn_t->ktime[i][j].c)
957 print_gputimes(fplog, k_log_str[i][j],
958 gpu_nbnxn_t->ktime[i][j].c, gpu_nbnxn_t->ktime[i][j].t, tot_gpu);
964 for (size_t k = 0; k < gtPME_EVENT_COUNT; k++)
966 if (gpu_pme_t->timing[k].c)
968 print_gputimes(fplog, PMEStageNames[k],
969 gpu_pme_t->timing[k].c,
970 gpu_pme_t->timing[k].t,
975 if (gpu_nbnxn_t->pruneTime.c)
977 print_gputimes(fplog, "Pruning kernel", gpu_nbnxn_t->pruneTime.c, gpu_nbnxn_t->pruneTime.t, tot_gpu);
979 print_gputimes(fplog, "F D2H", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_d2h_t, tot_gpu);
980 fprintf(fplog, "%s\n", hline);
981 print_gputimes(fplog, "Total ", gpu_nbnxn_t->nb_c, tot_gpu, tot_gpu);
982 fprintf(fplog, "%s\n", hline);
983 if (gpu_nbnxn_t->dynamicPruneTime.c)
985 /* We print the dynamic pruning kernel timings after a separator
986 * and avoid adding it to tot_gpu as this is not in the force
987 * overlap. We print the fraction as relative to the rest.
989 print_gputimes(fplog, "*Dynamic pruning", gpu_nbnxn_t->dynamicPruneTime.c, gpu_nbnxn_t->dynamicPruneTime.t, tot_gpu);
990 fprintf(fplog, "%s\n", hline);
992 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
993 if (gpu_nbnxn_t->nb_c > 0 && wc->wcc[ewcFORCE].n > 0)
995 fprintf(fplog, "\nAverage per-step force GPU/CPU evaluation time ratio: %.3f ms/%.3f ms = %.3f\n",
996 tot_gpu/gpu_nbnxn_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
1000 /* only print notes related to CPU-GPU load balance with PME */
1001 if (wc->wcc[ewcPMEMESH].n > 0)
1003 fprintf(fplog, "For optimal resource utilization this ratio should be close to 1\n");
1005 /* print note if the imbalance is high with PME case in which
1006 * CPU-GPU load balancing is possible */
1007 if (gpu_cpu_ratio < 0.8 || gpu_cpu_ratio > 1.25)
1009 /* Only the sim master calls this function, so always print to stderr */
1010 if (gpu_cpu_ratio < 0.8)
1014 /* The user could have used -notunepme,
1015 * but we currently can't check that here.
1017 GMX_LOG(mdlog.warning).asParagraph().appendText(
1018 "NOTE: The CPU has >25% more load than the GPU. This imbalance wastes\n"
1019 " GPU resources. Maybe the domain decomposition limits the PME tuning.\n"
1020 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
1024 /* We should not end up here, unless the box is
1025 * too small for increasing the cut-off for PME tuning.
1027 GMX_LOG(mdlog.warning).asParagraph().appendText(
1028 "NOTE: The CPU has >25% more load than the GPU. This imbalance wastes\n"
1032 if (gpu_cpu_ratio > 1.25)
1034 GMX_LOG(mdlog.warning).asParagraph().appendText(
1035 "NOTE: The GPU has >25% more load than the CPU. This imbalance wastes\n"
1044 GMX_LOG(mdlog.warning).asParagraph().appendText(
1045 "MPI_Barrier was called before each cycle start/stop\n"
1046 "call, so timings are not those of real runs.");
1049 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
1050 (cyc_sum[ewcDOMDEC] > tot*0.1 ||
1051 cyc_sum[ewcNS] > tot*0.1))
1053 /* Only the sim master calls this function, so always print to stderr */
1054 if (wc->wcc[ewcDOMDEC].n == 0)
1056 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1057 "NOTE: %d %% of the run time was spent in pair search,\n"
1058 " you might want to increase nstlist (this has no effect on accuracy)\n",
1059 gmx::roundToInt(100*cyc_sum[ewcNS]/tot));
1063 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1064 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
1065 " %d %% of the run time was spent in pair search,\n"
1066 " you might want to increase nstlist (this has no effect on accuracy)\n",
1067 gmx::roundToInt(100*cyc_sum[ewcDOMDEC]/tot),
1068 gmx::roundToInt(100*cyc_sum[ewcNS]/tot));
1072 if (cyc_sum[ewcMoveE] > tot*0.05)
1074 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1075 "NOTE: %d %% of the run time was spent communicating energies,\n"
1076 " you might want to increase some nst* mdp options\n",
1077 gmx::roundToInt(100*cyc_sum[ewcMoveE]/tot));
1081 extern int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
1088 return wc->reset_counters;
1091 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, int64_t reset_counters)
1098 wc->reset_counters = reset_counters;
1101 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
1103 if (useCycleSubcounters && wc != nullptr)
1105 wc->wcsc[ewcs].start = gmx_cycles_read();
1109 void wallcycle_sub_start_nocount(gmx_wallcycle_t wc, int ewcs)
1111 if (useCycleSubcounters && wc != nullptr)
1113 wallcycle_sub_start(wc, ewcs);
1118 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
1120 if (useCycleSubcounters && wc != nullptr)
1122 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;