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,2021, 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 //! Whether wallcycle debugging is enabled
62 constexpr bool gmx_unused enableWallcycleDebug = (DEBUG_WCYCLE != 0);
63 //! True if only the master rank should print debugging output
64 constexpr bool gmx_unused onlyMasterDebugPrints = true;
65 //! True if cycle counter nesting depth debuggin prints are enabled
66 constexpr bool gmx_unused debugPrintDepth = false /* enableWallcycleDebug */;
69 # include "gromacs/utility/fatalerror.h"
72 /* Each name should not exceed 19 printing characters
73 (ie. terminating null can be twentieth) */
74 static const char* wcn[ewcNR] = { "Run",
97 "Wait PME GPU spread",
99 "PME solve", /* the strings for FFT/solve are repeated here for mixed mode counters */
100 "Wait PME GPU gather",
103 "Wait GPU NB nonloc.",
105 "Wait GPU state copy",
106 "NB X/F buffer ops.",
120 static const char* wcsn[ewcsNR] = {
131 "NS search non-loc.",
135 "Listed buffer ops.",
137 "Nonbonded F kernel",
140 "Launch NB GPU tasks",
141 "Launch Bonded GPU tasks",
142 "Launch PME GPU tasks",
144 "Ewald F correction",
147 "Clear force buffer",
148 "Launch GPU NB X buffer ops.",
149 "Launch GPU NB F buffer ops.",
150 "Launch GPU Comm. coord.",
151 "Launch GPU Comm. force.",
156 /* PME GPU timing events' names - correspond to the enum in the gpu_timing.h */
157 static const char* enumValuetoString(PmeStage enumValue)
159 constexpr gmx::EnumerationArray<PmeStage, const char*> pmeStageNames = {
160 "PME spline", "PME spread", "PME spline + spread", "PME 3D-FFT r2c",
161 "PME solve", "PME 3D-FFT c2r", "PME gather"
163 return pmeStageNames[enumValue];
166 bool wallcycle_have_counter()
168 return gmx_cycles_have_counter();
171 gmx_wallcycle_t wallcycle_init(FILE* fplog, int resetstep, const t_commrec* cr)
176 if (!wallcycle_have_counter())
183 wc->haveInvalidCount = FALSE;
184 wc->wc_barrier = FALSE;
185 wc->wcc_all = nullptr;
188 wc->reset_counters = resetstep;
193 if (cr != nullptr && PAR(cr) && getenv("GMX_CYCLE_BARRIER") != nullptr)
197 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
199 wc->wc_barrier = TRUE;
203 snew(wc->wcc, ewcNR);
204 if (getenv("GMX_CYCLE_ALL") != nullptr)
208 fprintf(fplog, "\nWill time all the code during the run\n\n");
210 snew(wc->wcc_all, ewcNR * ewcNR);
213 if (sc_useCycleSubcounters)
215 snew(wc->wcsc, ewcsNR);
220 wc->isMasterRank = MASTER(cr);
226 void wallcycle_destroy(gmx_wallcycle_t wc)
233 if (wc->wcc != nullptr)
237 if (wc->wcc_all != nullptr)
241 if (wc->wcsc != nullptr)
249 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
251 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
253 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d", wc->count_depth + 1);
255 wc->counterlist[wc->count_depth] = ewc;
258 if (debugPrintDepth && (!onlyMasterDebugPrints || wc->isMasterRank))
260 std::string indentStr(4 * wc->count_depth, ' ');
261 fprintf(stderr, "%swcycle_start depth %d, %s\n", indentStr.c_str(), wc->count_depth, wcn[ewc]);
265 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
267 if (debugPrintDepth && (!onlyMasterDebugPrints || wc->isMasterRank))
269 std::string indentStr(4 * wc->count_depth, ' ');
270 fprintf(stderr, "%swcycle_stop depth %d, %s\n", indentStr.c_str(), wc->count_depth, wcn[ewc]);
275 if (wc->count_depth < 0)
277 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
279 if (wc->counterlist[wc->count_depth] != ewc)
282 "wallcycle mismatch at stop, start %s, stop %s",
283 wcn[wc->counterlist[wc->count_depth]],
289 void wallcycle_get(gmx_wallcycle_t wc, int ewc, int* n, double* c)
292 *c = static_cast<double>(wc->wcc[ewc].c);
295 void wallcycle_sub_get(gmx_wallcycle_t wc, int ewcs, int* n, double* c)
297 if (sc_useCycleSubcounters && wc != nullptr)
299 *n = wc->wcsc[ewcs].n;
300 *c = static_cast<double>(wc->wcsc[ewcs].c);
304 void wallcycle_reset_all(gmx_wallcycle_t wc)
313 for (i = 0; i < ewcNR; i++)
318 wc->haveInvalidCount = FALSE;
322 for (i = 0; i < ewcNR * ewcNR; i++)
324 wc->wcc_all[i].n = 0;
325 wc->wcc_all[i].c = 0;
330 for (i = 0; i < ewcsNR; i++)
338 static gmx_bool is_pme_counter(int ewc)
340 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
343 static gmx_bool is_pme_subcounter(int ewc)
345 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
348 void wallcycleBarrier(gmx_wallcycle* wc)
353 MPI_Barrier(wc->cr->mpi_comm_mygroup);
356 GMX_UNUSED_VALUE(wc);
360 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
361 static void subtract_cycles(wallcc_t* wcc, int ewc_main, int ewc_sub)
363 if (wcc[ewc_sub].n > 0)
365 if (wcc[ewc_main].c >= wcc[ewc_sub].c)
367 wcc[ewc_main].c -= wcc[ewc_sub].c;
371 /* Something is wrong with the cycle counting */
377 void wallcycle_scale_by_num_threads(gmx_wallcycle_t wc, bool isPmeRank, int nthreads_pp, int nthreads_pme)
384 for (int i = 0; i < ewcNR; i++)
386 if (is_pme_counter(i) || (i == ewcRUN && isPmeRank))
388 wc->wcc[i].c *= nthreads_pme;
392 for (int j = 0; j < ewcNR; j++)
394 wc->wcc_all[i * ewcNR + j].c *= nthreads_pme;
400 wc->wcc[i].c *= nthreads_pp;
404 for (int j = 0; j < ewcNR; j++)
406 wc->wcc_all[i * ewcNR + j].c *= nthreads_pp;
411 if (sc_useCycleSubcounters && wc->wcsc && !isPmeRank)
413 for (int i = 0; i < ewcsNR; i++)
415 wc->wcsc[i].c *= nthreads_pp;
420 /* TODO Make an object for this function to return, containing some
421 * vectors of something like wallcc_t for the summed wcc, wcc_all and
422 * wcsc, AND the original wcc for rank 0.
424 * The GPU timing is reported only for rank 0, so we want to preserve
425 * the original wcycle on that rank. Rank 0 also reports the global
426 * counts before that, so needs something to contain the global data
427 * without over-writing the rank-0 data. The current implementation
428 * uses cycles_sum to manage this, which works OK now because wcsc and
429 * wcc_all are unused by the GPU reporting, but it is not satisfactory
430 * for the future. Also, there's no need for MPI_Allreduce, since
431 * only MASTERRANK uses any of the results. */
432 WallcycleCounts wallcycle_sum(const t_commrec* cr, gmx_wallcycle_t wc)
434 WallcycleCounts cycles_sum;
436 double cycles[int(ewcNR) + int(ewcsNR)];
438 double cycles_n[int(ewcNR) + int(ewcsNR) + 1];
445 /* Default construction of std::array of non-class T can leave
446 the values indeterminate, just like a C array */
453 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
454 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
456 subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
458 if (cr->npmenodes == 0)
460 /* All nodes do PME (or no PME at all) */
461 subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
465 /* The are PME-only nodes */
466 if (wcc[ewcPMEMESH].n > 0)
468 /* This must be a PME only node, calculate the Wait + Comm. time */
469 GMX_ASSERT(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c,
470 "Total run ticks must be greater than PME-only ticks");
471 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
475 /* Store the cycles in a double buffer for summing */
476 for (i = 0; i < ewcNR; i++)
479 cycles_n[i] = static_cast<double>(wcc[i].n);
481 cycles[i] = static_cast<double>(wcc[i].c);
486 for (i = 0; i < ewcsNR; i++)
489 cycles_n[ewcNR + i] = static_cast<double>(wc->wcsc[i].n);
491 cycles[ewcNR + i] = static_cast<double>(wc->wcsc[i].c);
499 double buf[int(ewcNR) + int(ewcsNR) + 1];
501 // TODO this code is used only at the end of the run, so we
502 // can just do a simple reduce of haveInvalidCount in
503 // wallcycle_print, and avoid bugs
504 cycles_n[nsum] = (wc->haveInvalidCount ? 1 : 0);
505 // TODO Use MPI_Reduce
506 MPI_Allreduce(cycles_n, buf, nsum + 1, MPI_DOUBLE, MPI_MAX, cr->mpi_comm_mysim);
507 for (i = 0; i < ewcNR; i++)
509 wcc[i].n = gmx::roundToInt(buf[i]);
511 wc->haveInvalidCount = (buf[nsum] > 0);
514 for (i = 0; i < ewcsNR; i++)
516 wc->wcsc[i].n = gmx::roundToInt(buf[ewcNR + i]);
520 // TODO Use MPI_Reduce
521 MPI_Allreduce(cycles, cycles_sum.data(), nsum, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
523 if (wc->wcc_all != nullptr)
525 double *buf_all, *cyc_all;
527 snew(cyc_all, ewcNR * ewcNR);
528 snew(buf_all, ewcNR * ewcNR);
529 for (i = 0; i < ewcNR * ewcNR; i++)
531 cyc_all[i] = wc->wcc_all[i].c;
533 // TODO Use MPI_Reduce
534 MPI_Allreduce(cyc_all, buf_all, ewcNR * ewcNR, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
535 for (i = 0; i < ewcNR * ewcNR; i++)
537 wc->wcc_all[i].c = static_cast<gmx_cycles_t>(buf_all[i]);
546 for (i = 0; i < nsum; i++)
548 cycles_sum[i] = cycles[i];
556 print_cycles(FILE* fplog, double c2t, const char* name, int nnodes, int nthreads, int ncalls, double c_sum, double tot)
558 char nnodes_str[STRLEN];
559 char nthreads_str[STRLEN];
560 char ncalls_str[STRLEN];
562 double percentage = (tot > 0.) ? (100. * c_sum / tot) : 0.;
568 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
571 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
575 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
579 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
583 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
592 /* Convert the cycle count to wallclock time for this task */
596 " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
607 static void print_gputimes(FILE* fplog, const char* name, int n, double t, double tot_t)
614 snprintf(num, sizeof(num), "%10d", n);
615 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t / n);
620 sprintf(avg_perf, " ");
622 if (t != tot_t && tot_t > 0)
624 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n", name, num, t / 1000, avg_perf, 100 * t / tot_t);
628 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n", name, "", t / 1000, avg_perf, 100.0);
632 static void print_header(FILE* fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
634 int nrank_tot = nrank_pp + nrank_pme;
637 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
640 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
642 /* Don't report doing PP+PME, because we can't tell here if
643 * this is RF, etc. */
647 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
650 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
652 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
655 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
659 fprintf(fplog, "\n\n");
660 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
661 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
665 void wallcycle_print(FILE* fplog,
666 const gmx::MDLogger& mdlog,
673 const WallcycleCounts& cyc_sum,
674 const gmx_wallclock_gpu_nbnxn_t* gpu_nbnxn_t,
675 const gmx_wallclock_gpu_pme_t* gpu_pme_t)
677 double tot, tot_for_pp, tot_for_rest, tot_cpu_overlap, gpu_cpu_ratio;
678 double c2t, c2t_pp, c2t_pme = 0;
679 int i, j, npp, nth_tot;
682 "-----------------------------------------------------------------------------";
689 GMX_ASSERT(nth_pp > 0, "Number of particle-particle threads must be >0");
690 GMX_ASSERT(nth_pme > 0, "Number of PME threads must be >0");
691 GMX_ASSERT(nnodes > 0, "Number of nodes must be >0");
692 GMX_ASSERT(npme >= 0, "Number of PME nodes cannot be negative");
694 /* npme is the number of PME-only ranks used, and we always do PP work */
695 GMX_ASSERT(npp > 0, "Number of particle-particle nodes must be >0");
697 nth_tot = npp * nth_pp + npme * nth_pme;
699 /* When using PME-only nodes, the next line is valid for both
700 PP-only and PME-only nodes because they started ewcRUN at the
702 tot = cyc_sum[ewcRUN];
707 /* TODO This is heavy handed, but until someone reworks the
708 code so that it is provably robust with respect to
709 non-positive values for all possible timer and cycle
710 counters, there is less value gained from printing whatever
711 timing data might still be sensible for some non-Jenkins
712 run, than is lost from diagnosing Jenkins FP exceptions on
713 runs about whose execution time we don't care. */
714 GMX_LOG(mdlog.warning)
716 .appendTextFormatted(
717 "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a "
723 if (wc->haveInvalidCount)
725 GMX_LOG(mdlog.warning)
728 "NOTE: Detected invalid cycle counts, probably because threads moved "
729 "between CPU cores that do not have synchronized cycle counters. Will not "
730 "print the cycle accounting.");
735 /* Conversion factor from cycles to seconds */
736 c2t = realtime / tot;
737 c2t_pp = c2t * nth_tot / static_cast<double>(npp * nth_pp);
740 c2t_pme = c2t * nth_tot / static_cast<double>(npme * nth_pme);
747 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");
749 print_header(fplog, npp, nth_pp, npme, nth_pme);
751 fprintf(fplog, "%s\n", hline);
752 for (i = ewcPPDURINGPME + 1; i < ewcNR; i++)
754 if (is_pme_subcounter(i))
756 /* Do not count these at all */
758 else if (npme > 0 && is_pme_counter(i))
760 /* Print timing information for PME-only nodes, but add an
761 * asterisk so the reader of the table can know that the
762 * walltimes are not meant to add up. The asterisk still
763 * fits in the required maximum of 19 characters. */
765 snprintf(buffer, STRLEN, "%s *", wcn[i]);
766 print_cycles(fplog, c2t_pme, buffer, npme, nth_pme, wc->wcc[i].n, cyc_sum[i], tot);
770 /* Print timing information when it is for a PP or PP+PME
772 print_cycles(fplog, c2t_pp, wcn[i], npp, nth_pp, wc->wcc[i].n, cyc_sum[i], tot);
773 tot_for_pp += cyc_sum[i];
776 if (wc->wcc_all != nullptr)
778 for (i = 0; i < ewcNR; i++)
780 for (j = 0; j < ewcNR; j++)
782 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
788 wc->wcc_all[i * ewcNR + j].n,
789 wc->wcc_all[i * ewcNR + j].c,
794 tot_for_rest = tot * npp * nth_pp / static_cast<double>(nth_tot);
795 print_cycles(fplog, c2t_pp, "Rest", npp, nth_pp, -1, tot_for_rest - tot_for_pp, tot);
796 fprintf(fplog, "%s\n", hline);
797 print_cycles(fplog, c2t, "Total", npp, nth_pp, -1, tot, tot);
798 fprintf(fplog, "%s\n", hline);
803 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
804 " twice the total reported, but the cycle count total and %% are correct.\n"
809 if (wc->wcc[ewcPMEMESH].n > 0)
811 // A workaround to not print breakdown when no subcounters were recorded.
812 // TODO: figure out and record PME GPU counters (what to do with the waiting ones?)
813 std::vector<int> validPmeSubcounterIndices;
814 for (i = ewcPPDURINGPME + 1; i < ewcNR; i++)
816 if (is_pme_subcounter(i) && wc->wcc[i].n > 0)
818 validPmeSubcounterIndices.push_back(i);
822 if (!validPmeSubcounterIndices.empty())
824 fprintf(fplog, " Breakdown of PME mesh computation\n");
825 fprintf(fplog, "%s\n", hline);
826 for (auto i : validPmeSubcounterIndices)
829 npme > 0 ? c2t_pme : c2t_pp,
831 npme > 0 ? npme : npp,
837 fprintf(fplog, "%s\n", hline);
841 if (sc_useCycleSubcounters && wc->wcsc)
843 fprintf(fplog, " Breakdown of PP computation\n");
844 fprintf(fplog, "%s\n", hline);
845 for (i = 0; i < ewcsNR; i++)
847 print_cycles(fplog, c2t_pp, wcsn[i], npp, nth_pp, wc->wcsc[i].n, cyc_sum[ewcNR + i], tot);
849 fprintf(fplog, "%s\n", hline);
852 /* print GPU timing summary */
853 double tot_gpu = 0.0;
856 for (auto key : keysOf(gpu_pme_t->timing))
858 tot_gpu += gpu_pme_t->timing[key].t;
863 const char* k_log_str[2][2] = { { "Nonbonded F kernel", "Nonbonded F+ene k." },
864 { "Nonbonded F+prune k.", "Nonbonded F+ene+prune k." } };
865 tot_gpu += gpu_nbnxn_t->pl_h2d_t + gpu_nbnxn_t->nb_h2d_t + gpu_nbnxn_t->nb_d2h_t;
867 /* add up the kernel timings */
868 for (i = 0; i < 2; i++)
870 for (j = 0; j < 2; j++)
872 tot_gpu += gpu_nbnxn_t->ktime[i][j].t;
875 tot_gpu += gpu_nbnxn_t->pruneTime.t;
877 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
878 if (wc->wcc[ewcPMEMESH].n > 0)
880 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
882 tot_cpu_overlap *= realtime * 1000 / tot; /* convert s to ms */
884 fprintf(fplog, "\n GPU timings\n%s\n", hline);
886 " Computing: Count Wall t (s) ms/step %c\n",
888 fprintf(fplog, "%s\n", hline);
889 print_gputimes(fplog, "Pair list H2D", gpu_nbnxn_t->pl_h2d_c, gpu_nbnxn_t->pl_h2d_t, tot_gpu);
890 print_gputimes(fplog, "X / q H2D", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_h2d_t, tot_gpu);
892 for (i = 0; i < 2; i++)
894 for (j = 0; j < 2; j++)
896 if (gpu_nbnxn_t->ktime[i][j].c)
898 print_gputimes(fplog,
900 gpu_nbnxn_t->ktime[i][j].c,
901 gpu_nbnxn_t->ktime[i][j].t,
908 for (auto key : keysOf(gpu_pme_t->timing))
910 if (gpu_pme_t->timing[key].c)
912 print_gputimes(fplog,
913 enumValuetoString(key),
914 gpu_pme_t->timing[key].c,
915 gpu_pme_t->timing[key].t,
920 if (gpu_nbnxn_t->pruneTime.c)
922 print_gputimes(fplog, "Pruning kernel", gpu_nbnxn_t->pruneTime.c, gpu_nbnxn_t->pruneTime.t, tot_gpu);
924 print_gputimes(fplog, "F D2H", gpu_nbnxn_t->nb_c, gpu_nbnxn_t->nb_d2h_t, tot_gpu);
925 fprintf(fplog, "%s\n", hline);
926 print_gputimes(fplog, "Total ", gpu_nbnxn_t->nb_c, tot_gpu, tot_gpu);
927 fprintf(fplog, "%s\n", hline);
928 if (gpu_nbnxn_t->dynamicPruneTime.c)
930 /* We print the dynamic pruning kernel timings after a separator
931 * and avoid adding it to tot_gpu as this is not in the force
932 * overlap. We print the fraction as relative to the rest.
934 print_gputimes(fplog,
936 gpu_nbnxn_t->dynamicPruneTime.c,
937 gpu_nbnxn_t->dynamicPruneTime.t,
939 fprintf(fplog, "%s\n", hline);
941 gpu_cpu_ratio = tot_gpu / tot_cpu_overlap;
942 if (gpu_nbnxn_t->nb_c > 0 && wc->wcc[ewcFORCE].n > 0)
945 "\nAverage per-step force GPU/CPU evaluation time ratio: %.3f ms/%.3f ms = "
947 tot_gpu / gpu_nbnxn_t->nb_c,
948 tot_cpu_overlap / wc->wcc[ewcFORCE].n,
952 /* only print notes related to CPU-GPU load balance with PME */
953 if (wc->wcc[ewcPMEMESH].n > 0)
955 fprintf(fplog, "For optimal resource utilization this ratio should be close to 1\n");
957 /* print note if the imbalance is high with PME case in which
958 * CPU-GPU load balancing is possible */
959 if (gpu_cpu_ratio < 0.8 || gpu_cpu_ratio > 1.25)
961 /* Only the sim master calls this function, so always print to stderr */
962 if (gpu_cpu_ratio < 0.8)
966 /* The user could have used -notunepme,
967 * but we currently can't check that here.
969 GMX_LOG(mdlog.warning)
972 "NOTE: The CPU has >25% more load than the GPU. This "
974 " GPU resources. Maybe the domain decomposition "
975 "limits the PME tuning.\n"
976 " In that case, try setting the DD grid manually "
977 "(-dd) or lowering -dds.");
981 /* We should not end up here, unless the box is
982 * too small for increasing the cut-off for PME tuning.
984 GMX_LOG(mdlog.warning)
987 "NOTE: The CPU has >25% more load than the GPU. This "
992 if (gpu_cpu_ratio > 1.25)
994 GMX_LOG(mdlog.warning)
997 "NOTE: The GPU has >25% more load than the CPU. This imbalance "
1007 GMX_LOG(mdlog.warning)
1010 "MPI_Barrier was called before each cycle start/stop\n"
1011 "call, so timings are not those of real runs.");
1014 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 && (cyc_sum[ewcDOMDEC] > tot * 0.1 || cyc_sum[ewcNS] > tot * 0.1))
1016 /* Only the sim master calls this function, so always print to stderr */
1017 if (wc->wcc[ewcDOMDEC].n == 0)
1019 GMX_LOG(mdlog.warning)
1021 .appendTextFormatted(
1022 "NOTE: %d %% of the run time was spent in pair search,\n"
1023 " you might want to increase nstlist (this has no effect on "
1025 gmx::roundToInt(100 * cyc_sum[ewcNS] / tot));
1029 GMX_LOG(mdlog.warning)
1031 .appendTextFormatted(
1032 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
1033 " %d %% of the run time was spent in pair search,\n"
1034 " you might want to increase nstlist (this has no effect on "
1036 gmx::roundToInt(100 * cyc_sum[ewcDOMDEC] / tot),
1037 gmx::roundToInt(100 * cyc_sum[ewcNS] / tot));
1041 if (cyc_sum[ewcMoveE] > tot * 0.05)
1043 GMX_LOG(mdlog.warning)
1045 .appendTextFormatted(
1046 "NOTE: %d %% of the run time was spent communicating energies,\n"
1047 " you might want to increase some nst* mdp options\n",
1048 gmx::roundToInt(100 * cyc_sum[ewcMoveE] / tot));
1052 extern int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
1059 return wc->reset_counters;
1062 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, int64_t reset_counters)
1069 wc->reset_counters = reset_counters;