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, 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"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/timing/cyclecounter.h"
49 #include "gromacs/timing/gpu_timing.h"
50 #include "gromacs/timing/wallcyclereporting.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/gmxassert.h"
53 #include "gromacs/utility/gmxmpi.h"
54 #include "gromacs/utility/logger.h"
55 #include "gromacs/utility/smalloc.h"
56 #include "gromacs/utility/snprintf.h"
58 static const bool useCycleSubcounters = GMX_CYCLE_SUBCOUNTERS;
60 /* DEBUG_WCYCLE adds consistency checking for the counters.
61 * It checks if you stop a counter different from the last
62 * one that was opened and if you do nest too deep.
64 /* #define DEBUG_WCYCLE */
67 #include "gromacs/utility/fatalerror.h"
77 typedef struct gmx_wallcycle
80 /* did we detect one or more invalid cycle counts */
81 gmx_bool haveInvalidCount;
82 /* variables for testing/debugging */
88 int counterlist[DEPTH_MAX];
92 gmx_cycles_t cycle_prev;
93 gmx_int64_t reset_counters;
95 MPI_Comm mpi_comm_mygroup;
100 /* Each name should not exceed 19 printing characters
101 (ie. terminating null can be twentieth) */
102 static const char *wcn[ewcNR] =
104 "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
105 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
106 "Comm. coord.", "Born radii", "Force", "Wait + Comm. F", "PME mesh",
107 "PME redist. X/F", "PME spread/gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve LJ", "PME solve Elec",
108 "PME wait for PP", "Wait + Recv. PME F", "Wait GPU nonlocal", "Wait GPU local", "NB X/F buffer ops.",
109 "Vsite spread", "COM pull force",
110 "Write traj.", "Update", "Constraints", "Comm. energies",
111 "Enforced rotation", "Add rot. forces", "Position swapping", "IMD", "Test"
114 static const char *wcsn[ewcsNR] =
116 "DD redist.", "DD NS grid + sort", "DD setup comm.",
117 "DD make top.", "DD make constr.", "DD top. other",
118 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
122 "Listed buffer ops.",
125 "Ewald F correction",
130 gmx_bool wallcycle_have_counter(void)
132 return gmx_cycles_have_counter();
135 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec gmx_unused *cr)
140 if (!wallcycle_have_counter())
147 wc->haveInvalidCount = FALSE;
148 wc->wc_barrier = FALSE;
149 wc->wcc_all = nullptr;
152 wc->reset_counters = resetstep;
155 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != nullptr)
159 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
161 wc->wc_barrier = TRUE;
162 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
166 snew(wc->wcc, ewcNR);
167 if (getenv("GMX_CYCLE_ALL") != nullptr)
171 fprintf(fplog, "\nWill time all the code during the run\n\n");
173 snew(wc->wcc_all, ewcNR*ewcNR);
176 if (useCycleSubcounters)
178 snew(wc->wcsc, ewcsNR);
188 void wallcycle_destroy(gmx_wallcycle_t wc)
195 if (wc->wcc != nullptr)
199 if (wc->wcc_all != nullptr)
203 if (wc->wcsc != nullptr)
210 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
213 wc->cycle_prev = cycle;
216 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
218 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
219 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
224 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
226 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
228 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
230 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
233 wc->counterlist[wc->count_depth] = ewc;
237 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
241 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
243 if (wc->count_depth < 0)
245 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
247 if (wc->counterlist[wc->count_depth] != ewc)
249 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
250 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
255 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
267 MPI_Barrier(wc->mpi_comm_mygroup);
272 debug_start_check(wc, ewc);
275 cycle = gmx_cycles_read();
276 wc->wcc[ewc].start = cycle;
277 if (wc->wcc_all != nullptr)
282 wallcycle_all_start(wc, ewc, cycle);
284 else if (wc->wc_depth == 3)
286 wallcycle_all_stop(wc, ewc, cycle);
291 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
298 wallcycle_start(wc, ewc);
302 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
304 gmx_cycles_t cycle, last;
314 MPI_Barrier(wc->mpi_comm_mygroup);
319 debug_stop_check(wc, ewc);
322 /* When processes or threads migrate between cores, the cycle counting
323 * can get messed up if the cycle counter on different cores are not
324 * synchronized. When this happens we expect both large negative and
325 * positive cycle differences. We can detect negative cycle differences.
326 * Detecting too large positive counts if difficult, since count can be
327 * large, especially for ewcRUN. If we detect a negative count,
328 * we will not print the cycle accounting table.
330 cycle = gmx_cycles_read();
331 if (cycle >= wc->wcc[ewc].start)
333 last = cycle - wc->wcc[ewc].start;
338 wc->haveInvalidCount = TRUE;
340 wc->wcc[ewc].c += last;
347 wallcycle_all_stop(wc, ewc, cycle);
349 else if (wc->wc_depth == 2)
351 wallcycle_all_start(wc, ewc, cycle);
358 void wallcycle_get(gmx_wallcycle_t wc, int ewc, int *n, double *c)
361 *c = static_cast<double>(wc->wcc[ewc].c);
364 void wallcycle_reset_all(gmx_wallcycle_t wc)
373 for (i = 0; i < ewcNR; i++)
378 wc->haveInvalidCount = FALSE;
382 for (i = 0; i < ewcNR*ewcNR; i++)
384 wc->wcc_all[i].n = 0;
385 wc->wcc_all[i].c = 0;
390 for (i = 0; i < ewcsNR; i++)
398 static gmx_bool is_pme_counter(int ewc)
400 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
403 static gmx_bool is_pme_subcounter(int ewc)
405 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
408 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
409 static void subtract_cycles(wallcc_t *wcc, int ewc_main, int ewc_sub)
411 if (wcc[ewc_sub].n > 0)
413 if (wcc[ewc_main].c >= wcc[ewc_sub].c)
415 wcc[ewc_main].c -= wcc[ewc_sub].c;
419 /* Something is wrong with the cycle counting */
425 void wallcycle_scale_by_num_threads(gmx_wallcycle_t wc, bool isPmeRank, int nthreads_pp, int nthreads_pme)
432 for (int i = 0; i < ewcNR; i++)
434 if (is_pme_counter(i) || (i == ewcRUN && isPmeRank))
436 wc->wcc[i].c *= nthreads_pme;
440 for (int j = 0; j < ewcNR; j++)
442 wc->wcc_all[i*ewcNR+j].c *= nthreads_pme;
448 wc->wcc[i].c *= nthreads_pp;
452 for (int j = 0; j < ewcNR; j++)
454 wc->wcc_all[i*ewcNR+j].c *= nthreads_pp;
459 if (useCycleSubcounters && wc->wcsc && !isPmeRank)
461 for (int i = 0; i < ewcsNR; i++)
463 wc->wcsc[i].c *= nthreads_pp;
468 /* TODO Make an object for this function to return, containing some
469 * vectors of something like wallcc_t for the summed wcc, wcc_all and
470 * wcsc, AND the original wcc for rank 0.
472 * The GPU timing is reported only for rank 0, so we want to preserve
473 * the original wcycle on that rank. Rank 0 also reports the global
474 * counts before that, so needs something to contain the global data
475 * without over-writing the rank-0 data. The current implementation
476 * uses cycles_sum to manage this, which works OK now because wcsc and
477 * wcc_all are unused by the GPU reporting, but it is not satisfactory
478 * for the future. Also, there's no need for MPI_Allreduce, since
479 * only MASTERRANK uses any of the results. */
480 WallcycleCounts wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
482 WallcycleCounts cycles_sum;
484 double cycles[ewcNR+ewcsNR];
486 double cycles_n[ewcNR+ewcsNR+1];
493 /* Default construction of std::array of non-class T can leave
494 the values indeterminate, just like a C array, and icc
502 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
503 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
505 subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
507 if (cr->npmenodes == 0)
509 /* All nodes do PME (or no PME at all) */
510 subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
514 /* The are PME-only nodes */
515 if (wcc[ewcPMEMESH].n > 0)
517 /* This must be a PME only node, calculate the Wait + Comm. time */
518 GMX_ASSERT(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c, "Total run ticks must be greater than PME-only ticks");
519 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
523 /* Store the cycles in a double buffer for summing */
524 for (i = 0; i < ewcNR; i++)
527 cycles_n[i] = static_cast<double>(wcc[i].n);
529 cycles[i] = static_cast<double>(wcc[i].c);
534 for (i = 0; i < ewcsNR; i++)
537 cycles_n[ewcNR+i] = static_cast<double>(wc->wcsc[i].n);
539 cycles[ewcNR+i] = static_cast<double>(wc->wcsc[i].c);
547 double buf[ewcNR+ewcsNR+1];
549 // TODO this code is used only at the end of the run, so we
550 // can just do a simple reduce of haveInvalidCount in
551 // wallcycle_print, and avoid bugs
552 cycles_n[nsum] = (wc->haveInvalidCount > 0 ? 1 : 0);
553 // TODO Use MPI_Reduce
554 MPI_Allreduce(cycles_n, buf, nsum + 1, MPI_DOUBLE, MPI_MAX,
556 for (i = 0; i < ewcNR; i++)
558 wcc[i].n = static_cast<int>(buf[i] + 0.5);
560 wc->haveInvalidCount = (buf[nsum] > 0);
563 for (i = 0; i < ewcsNR; i++)
565 wc->wcsc[i].n = static_cast<int>(buf[ewcNR+i] + 0.5);
569 // TODO Use MPI_Reduce
570 MPI_Allreduce(cycles, cycles_sum.data(), nsum, MPI_DOUBLE, MPI_SUM,
573 if (wc->wcc_all != nullptr)
575 double *buf_all, *cyc_all;
577 snew(cyc_all, ewcNR*ewcNR);
578 snew(buf_all, ewcNR*ewcNR);
579 for (i = 0; i < ewcNR*ewcNR; i++)
581 cyc_all[i] = wc->wcc_all[i].c;
583 // TODO Use MPI_Reduce
584 MPI_Allreduce(cyc_all, buf_all, ewcNR*ewcNR, MPI_DOUBLE, MPI_SUM,
586 for (i = 0; i < ewcNR*ewcNR; i++)
588 wc->wcc_all[i].c = static_cast<gmx_cycles_t>(buf_all[i]);
597 for (i = 0; i < nsum; i++)
599 cycles_sum[i] = cycles[i];
606 static void print_cycles(FILE *fplog, double c2t, const char *name,
607 int nnodes, int nthreads,
608 int ncalls, double c_sum, double tot)
610 char nnodes_str[STRLEN];
611 char nthreads_str[STRLEN];
612 char ncalls_str[STRLEN];
614 double percentage = (tot > 0.) ? (100. * c_sum / tot) : 0.;
620 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
623 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
627 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
631 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
635 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
644 /* Convert the cycle count to wallclock time for this task */
647 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
648 name, nnodes_str, nthreads_str, ncalls_str, wallt,
649 c_sum*1e-9, percentage);
653 static void print_gputimes(FILE *fplog, const char *name,
654 int n, double t, double tot_t)
661 snprintf(num, sizeof(num), "%10d", n);
662 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
667 sprintf(avg_perf, " ");
669 if (t != tot_t && tot_t > 0)
671 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
672 name, num, t/1000, avg_perf, 100 * t/tot_t);
676 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
677 name, "", t/1000, avg_perf, 100.0);
681 static void print_header(FILE *fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
683 int nrank_tot = nrank_pp + nrank_pme;
686 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
689 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
691 /* Don't report doing PP+PME, because we can't tell here if
692 * this is RF, etc. */
696 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
699 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
701 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
704 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
708 fprintf(fplog, "\n\n");
709 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
710 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
714 void wallcycle_print(FILE *fplog, const gmx::MDLogger &mdlog, int nnodes, int npme,
715 int nth_pp, int nth_pme, double realtime,
716 gmx_wallcycle_t wc, const WallcycleCounts &cyc_sum,
717 struct gmx_wallclock_gpu_t *gpu_t)
719 double tot, tot_for_pp, tot_for_rest, tot_gpu, tot_cpu_overlap, gpu_cpu_ratio, tot_k;
720 double c2t, c2t_pp, c2t_pme = 0;
721 int i, j, npp, nth_tot;
723 const char *hline = "-----------------------------------------------------------------------------";
730 GMX_ASSERT(nth_pp > 0, "Number of particle-particle threads must be >0");
731 GMX_ASSERT(nth_pme > 0, "Number of PME threads must be >0");
732 GMX_ASSERT(nnodes > 0, "Number of nodes must be >0");
733 GMX_ASSERT(npme >= 0, "Number of PME nodes cannot be negative");
735 /* npme is the number of PME-only ranks used, and we always do PP work */
736 GMX_ASSERT(npp > 0, "Number of particle-particle nodes must be >0");
738 nth_tot = npp*nth_pp + npme*nth_pme;
740 /* When using PME-only nodes, the next line is valid for both
741 PP-only and PME-only nodes because they started ewcRUN at the
743 tot = cyc_sum[ewcRUN];
748 /* TODO This is heavy handed, but until someone reworks the
749 code so that it is provably robust with respect to
750 non-positive values for all possible timer and cycle
751 counters, there is less value gained from printing whatever
752 timing data might still be sensible for some non-Jenkins
753 run, than is lost from diagnosing Jenkins FP exceptions on
754 runs about whose execution time we don't care. */
755 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
756 "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a time accounting",
761 if (wc->haveInvalidCount)
763 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.");
768 /* Conversion factor from cycles to seconds */
770 c2t_pp = c2t * nth_tot / static_cast<double>(npp*nth_pp);
773 c2t_pme = c2t * nth_tot / static_cast<double>(npme*nth_pme);
780 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");
782 print_header(fplog, npp, nth_pp, npme, nth_pme);
784 fprintf(fplog, "%s\n", hline);
785 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
787 if (is_pme_subcounter(i))
789 /* Do not count these at all */
791 else if (npme > 0 && is_pme_counter(i))
793 /* Print timing information for PME-only nodes, but add an
794 * asterisk so the reader of the table can know that the
795 * walltimes are not meant to add up. The asterisk still
796 * fits in the required maximum of 19 characters. */
798 snprintf(buffer, STRLEN, "%s *", wcn[i]);
799 print_cycles(fplog, c2t_pme, buffer,
801 wc->wcc[i].n, cyc_sum[i], tot);
805 /* Print timing information when it is for a PP or PP+PME
807 print_cycles(fplog, c2t_pp, wcn[i],
809 wc->wcc[i].n, cyc_sum[i], tot);
810 tot_for_pp += cyc_sum[i];
813 if (wc->wcc_all != nullptr)
815 for (i = 0; i < ewcNR; i++)
817 for (j = 0; j < ewcNR; j++)
819 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
820 print_cycles(fplog, c2t_pp, buf,
822 wc->wcc_all[i*ewcNR+j].n,
823 wc->wcc_all[i*ewcNR+j].c,
828 tot_for_rest = tot * npp * nth_pp / static_cast<double>(nth_tot);
829 print_cycles(fplog, c2t_pp, "Rest",
831 -1, tot_for_rest - tot_for_pp, tot);
832 fprintf(fplog, "%s\n", hline);
833 print_cycles(fplog, c2t, "Total",
836 fprintf(fplog, "%s\n", hline);
841 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
842 " twice the total reported, but the cycle count total and %% are correct.\n"
846 if (wc->wcc[ewcPMEMESH].n > 0)
848 fprintf(fplog, " Breakdown of PME mesh computation\n");
849 fprintf(fplog, "%s\n", hline);
850 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
852 if (is_pme_subcounter(i))
854 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i],
855 npme > 0 ? npme : npp, nth_pme,
856 wc->wcc[i].n, cyc_sum[i], tot);
859 fprintf(fplog, "%s\n", hline);
862 if (useCycleSubcounters && wc->wcsc)
864 fprintf(fplog, " Breakdown of PP computation\n");
865 fprintf(fplog, "%s\n", hline);
866 for (i = 0; i < ewcsNR; i++)
868 print_cycles(fplog, c2t_pp, wcsn[i],
870 wc->wcsc[i].n, cyc_sum[ewcNR+i], tot);
872 fprintf(fplog, "%s\n", hline);
875 /* print GPU timing summary */
878 const char *k_log_str[2][2] = {
879 {"Nonbonded F kernel", "Nonbonded F+ene k."},
880 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
883 tot_gpu = gpu_t->pl_h2d_t + gpu_t->nb_h2d_t + gpu_t->nb_d2h_t;
885 /* add up the kernel timings */
887 for (i = 0; i < 2; i++)
889 for (j = 0; j < 2; j++)
891 tot_k += gpu_t->ktime[i][j].t;
894 tot_gpu += tot_k + gpu_t->pruneTime.t;
896 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
897 if (wc->wcc[ewcPMEMESH].n > 0)
899 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
901 tot_cpu_overlap *= realtime*1000/tot; /* convert s to ms */
903 fprintf(fplog, "\n GPU timings\n%s\n", hline);
904 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
905 fprintf(fplog, "%s\n", hline);
906 print_gputimes(fplog, "Pair list H2D",
907 gpu_t->pl_h2d_c, gpu_t->pl_h2d_t, tot_gpu);
908 print_gputimes(fplog, "X / q H2D",
909 gpu_t->nb_c, gpu_t->nb_h2d_t, tot_gpu);
911 for (i = 0; i < 2; i++)
913 for (j = 0; j < 2; j++)
915 if (gpu_t->ktime[i][j].c)
917 print_gputimes(fplog, k_log_str[i][j],
918 gpu_t->ktime[i][j].c, gpu_t->ktime[i][j].t, tot_gpu);
922 if (gpu_t->pruneTime.c)
924 print_gputimes(fplog, "Pruning kernel", gpu_t->pruneTime.c, gpu_t->pruneTime.t, tot_gpu);
926 print_gputimes(fplog, "F D2H", gpu_t->nb_c, gpu_t->nb_d2h_t, tot_gpu);
927 fprintf(fplog, "%s\n", hline);
928 print_gputimes(fplog, "Total ", gpu_t->nb_c, tot_gpu, tot_gpu);
929 fprintf(fplog, "%s\n", hline);
930 if (gpu_t->dynamicPruneTime.c)
932 /* We print the dynamic pruning kernel timings after a separator
933 * and avoid adding it to tot_gpu as this is not in the force
934 * overlap. We print the fraction as relative to the rest.
936 print_gputimes(fplog, "*Dynamic pruning", gpu_t->dynamicPruneTime.c, gpu_t->dynamicPruneTime.t, tot_gpu);
937 fprintf(fplog, "%s\n", hline);
940 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
941 if (gpu_t->nb_c > 0 && wc->wcc[ewcFORCE].n > 0)
943 fprintf(fplog, "\nAverage per-step force GPU/CPU evaluation time ratio: %.3f ms/%.3f ms = %.3f\n",
944 tot_gpu/gpu_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
948 /* only print notes related to CPU-GPU load balance with PME */
949 if (wc->wcc[ewcPMEMESH].n > 0)
951 fprintf(fplog, "For optimal performance this ratio should be close to 1!\n");
953 /* print note if the imbalance is high with PME case in which
954 * CPU-GPU load balancing is possible */
955 if (gpu_cpu_ratio < 0.75 || gpu_cpu_ratio > 1.2)
957 /* Only the sim master calls this function, so always print to stderr */
958 if (gpu_cpu_ratio < 0.75)
962 /* The user could have used -notunepme,
963 * but we currently can't check that here.
965 GMX_LOG(mdlog.warning).asParagraph().appendText(
966 "NOTE: The GPU has >25% less load than the CPU. This imbalance causes\n"
967 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
968 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
972 /* We should not end up here, unless the box is
973 * too small for increasing the cut-off for PME tuning.
975 GMX_LOG(mdlog.warning).asParagraph().appendText(
976 "NOTE: The GPU has >25% less load than the CPU. This imbalance causes\n"
977 " performance loss.");
980 if (gpu_cpu_ratio > 1.2)
982 GMX_LOG(mdlog.warning).asParagraph().appendText(
983 "NOTE: The GPU has >20% more load than the CPU. This imbalance causes\n"
984 " performance loss, consider using a shorter cut-off and a finer PME grid.");
992 GMX_LOG(mdlog.warning).asParagraph().appendText(
993 "MPI_Barrier was called before each cycle start/stop\n"
994 "call, so timings are not those of real runs.");
997 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
998 (cyc_sum[ewcDOMDEC] > tot*0.1 ||
999 cyc_sum[ewcNS] > tot*0.1))
1001 /* Only the sim master calls this function, so always print to stderr */
1002 if (wc->wcc[ewcDOMDEC].n == 0)
1004 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1005 "NOTE: %d %% of the run time was spent in pair search,\n"
1006 " you might want to increase nstlist (this has no effect on accuracy)\n",
1007 (int)(100*cyc_sum[ewcNS]/tot+0.5));
1011 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1012 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
1013 " %d %% of the run time was spent in pair search,\n"
1014 " you might want to increase nstlist (this has no effect on accuracy)\n",
1015 (int)(100*cyc_sum[ewcDOMDEC]/tot+0.5),
1016 (int)(100*cyc_sum[ewcNS]/tot+0.5));
1020 if (cyc_sum[ewcMoveE] > tot*0.05)
1022 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1023 "NOTE: %d %% of the run time was spent communicating energies,\n"
1024 " you might want to use the -gcom option of mdrun\n",
1025 (int)(100*cyc_sum[ewcMoveE]/tot+0.5));
1029 extern gmx_int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
1036 return wc->reset_counters;
1039 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_int64_t reset_counters)
1046 wc->reset_counters = reset_counters;
1049 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
1051 if (useCycleSubcounters && wc != nullptr)
1053 wc->wcsc[ewcs].start = gmx_cycles_read();
1057 void wallcycle_sub_start_nocount(gmx_wallcycle_t wc, int ewcs)
1059 if (useCycleSubcounters && wc != nullptr)
1061 wallcycle_sub_start(wc, ewcs);
1066 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
1068 if (useCycleSubcounters && wc != nullptr)
1070 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;