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, 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"
45 #include "gromacs/legacyheaders/md_logging.h"
46 #include "gromacs/legacyheaders/types/commrec.h"
47 #include "gromacs/timing/cyclecounter.h"
48 #include "gromacs/timing/gpu_timing.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/gmxmpi.h"
52 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/utility/snprintf.h"
55 /* DEBUG_WCYCLE adds consistency checking for the counters.
56 * It checks if you stop a counter different from the last
57 * one that was opened and if you do nest too deep.
59 /* #define DEBUG_WCYCLE */
62 #include "gromacs/utility/fatalerror.h"
73 typedef struct gmx_wallcycle
76 /* variables for testing/debugging */
82 int counterlist[DEPTH_MAX];
86 gmx_cycles_t cycle_prev;
87 gmx_int64_t reset_counters;
89 MPI_Comm mpi_comm_mygroup;
93 #ifdef GMX_CYCLE_SUBCOUNTERS
99 /* Each name should not exceed 19 printing characters
100 (ie. terminating null can be twentieth) */
101 static const char *wcn[ewcNR] =
103 "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
104 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
105 "Comm. coord.", "Born radii", "Force", "Wait + Comm. F", "PME mesh",
106 "PME redist. X/F", "PME spread/gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve LJ", "PME solve Elec",
107 "PME wait for PP", "Wait + Recv. PME F", "Wait GPU nonlocal", "Wait GPU local", "Wait GPU loc. est.", "NB X/F buffer ops.",
108 "Vsite spread", "COM pull force",
109 "Write traj.", "Update", "Constraints", "Comm. energies",
110 "Enforced rotation", "Add rot. forces", "Coordinate swapping", "IMD", "Test"
113 #ifdef GMX_CYCLE_SUBCOUNTERS
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.",
124 "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,
136 int nthreads_pp, int nthreads_pme)
141 if (!wallcycle_have_counter())
148 wc->wc_barrier = FALSE;
152 wc->reset_counters = resetstep;
153 wc->nthreads_pp = nthreads_pp;
154 wc->nthreads_pme = nthreads_pme;
155 wc->cycles_sum = NULL;
158 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != NULL)
162 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
164 wc->wc_barrier = TRUE;
165 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
169 snew(wc->wcc, ewcNR);
170 if (getenv("GMX_CYCLE_ALL") != NULL)
174 fprintf(fplog, "\nWill time all the code during the run\n\n");
176 snew(wc->wcc_all, ewcNR*ewcNR);
179 #ifdef GMX_CYCLE_SUBCOUNTERS
180 snew(wc->wcsc, ewcsNR);
190 void wallcycle_destroy(gmx_wallcycle_t wc)
201 if (wc->wcc_all != NULL)
205 #ifdef GMX_CYCLE_SUBCOUNTERS
206 if (wc->wcsc != NULL)
214 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
217 wc->cycle_prev = cycle;
220 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
222 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
223 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
228 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
230 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
232 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
234 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
237 wc->counterlist[wc->count_depth] = ewc;
241 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
245 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
247 if (wc->count_depth < 0)
249 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
251 if (wc->counterlist[wc->count_depth] != ewc)
253 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
254 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
259 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
271 MPI_Barrier(wc->mpi_comm_mygroup);
276 debug_start_check(wc, ewc);
279 cycle = gmx_cycles_read();
280 wc->wcc[ewc].start = cycle;
281 if (wc->wcc_all != NULL)
286 wallcycle_all_start(wc, ewc, cycle);
288 else if (wc->wc_depth == 3)
290 wallcycle_all_stop(wc, ewc, cycle);
295 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
302 wallcycle_start(wc, ewc);
306 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
308 gmx_cycles_t cycle, last;
318 MPI_Barrier(wc->mpi_comm_mygroup);
323 debug_stop_check(wc, ewc);
326 cycle = gmx_cycles_read();
327 last = cycle - wc->wcc[ewc].start;
328 wc->wcc[ewc].c += last;
335 wallcycle_all_stop(wc, ewc, cycle);
337 else if (wc->wc_depth == 2)
339 wallcycle_all_start(wc, ewc, cycle);
346 void wallcycle_get(gmx_wallcycle_t wc, int ewc, int *n, double *c)
349 *c = static_cast<double>(wc->wcc[ewc].c);
352 void wallcycle_reset_all(gmx_wallcycle_t wc)
361 for (i = 0; i < ewcNR; i++)
368 for (i = 0; i < ewcNR*ewcNR; i++)
370 wc->wcc_all[i].n = 0;
371 wc->wcc_all[i].c = 0;
374 #ifdef GMX_CYCLE_SUBCOUNTERS
375 for (i = 0; i < ewcsNR; i++)
383 static gmx_bool is_pme_counter(int ewc)
385 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
388 static gmx_bool is_pme_subcounter(int ewc)
390 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
393 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
394 static void subtract_cycles(wallcc_t *wcc, int ewc_main, int ewc_sub)
396 if (wcc[ewc_sub].n > 0)
398 GMX_ASSERT(wcc[ewc_main].c >= wcc[ewc_sub].c, "Subcounter cannot have more ticks than parent");
400 wcc[ewc_main].c -= wcc[ewc_sub].c;
404 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
407 double cycles[ewcNR+ewcsNR];
409 double cycles_n[ewcNR+ewcsNR];
410 double buf[ewcNR+ewcsNR];
411 double *buf_all, *cyc_all;
421 snew(wc->cycles_sum, ewcNR+ewcsNR);
425 /* The GPU wait estimate counter is used for load balancing only
426 * and will mess up the total due to double counting: clear it.
428 wcc[ewcWAIT_GPU_NB_L_EST].n = 0;
429 wcc[ewcWAIT_GPU_NB_L_EST].c = 0;
431 for (i = 0; i < ewcNR; i++)
433 if (is_pme_counter(i) || (i == ewcRUN && cr->duty == DUTY_PME))
435 wcc[i].c *= wc->nthreads_pme;
439 for (j = 0; j < ewcNR; j++)
441 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pme;
447 wcc[i].c *= wc->nthreads_pp;
451 for (j = 0; j < ewcNR; j++)
453 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pp;
459 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
460 subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
462 subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
464 if (cr->npmenodes == 0)
466 /* All nodes do PME (or no PME at all) */
467 subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
471 /* The are PME-only nodes */
472 if (wcc[ewcPMEMESH].n > 0)
474 /* This must be a PME only node, calculate the Wait + Comm. time */
475 GMX_ASSERT(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c, "Total run ticks must be greater than PME-only ticks");
476 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
480 /* Store the cycles in a double buffer for summing */
481 for (i = 0; i < ewcNR; i++)
484 cycles_n[i] = static_cast<double>(wcc[i].n);
486 cycles[i] = static_cast<double>(wcc[i].c);
489 #ifdef GMX_CYCLE_SUBCOUNTERS
490 for (i = 0; i < ewcsNR; i++)
492 wc->wcsc[i].c *= wc->nthreads_pp;
494 cycles_n[ewcNR+i] = static_cast<double>(wc->wcsc[i].n);
496 cycles[ewcNR+i] = static_cast<double>(wc->wcsc[i].c);
504 MPI_Allreduce(cycles_n, buf, nsum, MPI_DOUBLE, MPI_MAX,
506 for (i = 0; i < ewcNR; i++)
508 wcc[i].n = static_cast<int>(buf[i] + 0.5);
510 #ifdef GMX_CYCLE_SUBCOUNTERS
511 for (i = 0; i < ewcsNR; i++)
513 wc->wcsc[i].n = static_cast<int>(buf[ewcNR+i] + 0.5);
517 MPI_Allreduce(cycles, wc->cycles_sum, nsum, MPI_DOUBLE, MPI_SUM,
520 if (wc->wcc_all != NULL)
522 snew(cyc_all, ewcNR*ewcNR);
523 snew(buf_all, ewcNR*ewcNR);
524 for (i = 0; i < ewcNR*ewcNR; i++)
526 cyc_all[i] = wc->wcc_all[i].c;
528 MPI_Allreduce(cyc_all, buf_all, ewcNR*ewcNR, MPI_DOUBLE, MPI_SUM,
530 for (i = 0; i < ewcNR*ewcNR; i++)
532 wc->wcc_all[i].c = static_cast<gmx_cycles_t>(buf_all[i]);
541 for (i = 0; i < nsum; i++)
543 wc->cycles_sum[i] = cycles[i];
548 static void print_cycles(FILE *fplog, double c2t, const char *name,
549 int nnodes, int nthreads,
550 int ncalls, double c_sum, double tot)
553 char nthreads_str[6];
556 double percentage = (tot > 0.) ? (100. * c_sum / tot) : 0.;
562 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
565 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
569 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
573 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
577 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
586 /* Convert the cycle count to wallclock time for this task */
589 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
590 name, nnodes_str, nthreads_str, ncalls_str, wallt,
591 c_sum*1e-9, percentage);
595 static void print_gputimes(FILE *fplog, const char *name,
596 int n, double t, double tot_t)
603 snprintf(num, sizeof(num), "%10d", n);
604 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
609 sprintf(avg_perf, " ");
611 if (t != tot_t && tot_t > 0)
613 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
614 name, num, t/1000, avg_perf, 100 * t/tot_t);
618 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
619 name, "", t/1000, avg_perf, 100.0);
623 static void print_header(FILE *fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
625 int nrank_tot = nrank_pp + nrank_pme;
628 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
631 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
633 /* Don't report doing PP+PME, because we can't tell here if
634 * this is RF, etc. */
638 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
641 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
643 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
646 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
650 fprintf(fplog, "\n\n");
651 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
652 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
655 void wallcycle_print(FILE *fplog, int nnodes, int npme, double realtime,
656 gmx_wallcycle_t wc, struct gmx_wallclock_gpu_t *gpu_t)
659 double tot, tot_for_pp, tot_for_rest, tot_gpu, tot_cpu_overlap, gpu_cpu_ratio, tot_k;
660 double c2t, c2t_pp, c2t_pme = 0;
661 int i, j, npp, nth_pp, nth_pme, nth_tot;
663 const char *hline = "-----------------------------------------------------------------------------";
670 nth_pp = wc->nthreads_pp;
671 GMX_ASSERT(nth_pp > 0, "Number of particle-particle threads must be >0");
673 nth_pme = wc->nthreads_pme;
674 GMX_ASSERT(nth_pme > 0, "Number of PME threads must be >0");
676 cyc_sum = wc->cycles_sum;
678 GMX_ASSERT(nnodes > 0, "Number of nodes must be >0");
679 GMX_ASSERT(npme >= 0, "Number of PME nodes cannot be negative");
681 /* npme is the number of PME-only ranks used, and we always do PP work */
682 GMX_ASSERT(npp > 0, "Number of particle-particle nodes must be >0");
684 nth_tot = npp*nth_pp + npme*nth_pme;
686 /* When using PME-only nodes, the next line is valid for both
687 PP-only and PME-only nodes because they started ewcRUN at the
689 tot = cyc_sum[ewcRUN];
694 /* TODO This is heavy handed, but until someone reworks the
695 code so that it is provably robust with respect to
696 non-positive values for all possible timer and cycle
697 counters, there is less value gained from printing whatever
698 timing data might still be sensible for some non-Jenkins
699 run, than is lost from diagnosing Jenkins FP exceptions on
700 runs about whose execution time we don't care. */
701 md_print_warn(NULL, fplog, "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a time accounting\n", tot);
705 /* Conversion factor from cycles to seconds */
707 c2t_pp = c2t * nth_tot / static_cast<double>(npp*nth_pp);
710 c2t_pme = c2t * nth_tot / static_cast<double>(npme*nth_pme);
717 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");
719 print_header(fplog, npp, nth_pp, npme, nth_pme);
721 fprintf(fplog, "%s\n", hline);
722 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
724 if (is_pme_subcounter(i))
726 /* Do not count these at all */
728 else if (npme > 0 && is_pme_counter(i))
730 /* Print timing information for PME-only nodes, but add an
731 * asterisk so the reader of the table can know that the
732 * walltimes are not meant to add up. The asterisk still
733 * fits in the required maximum of 19 characters. */
735 snprintf(buffer, STRLEN, "%s *", wcn[i]);
736 print_cycles(fplog, c2t_pme, buffer,
738 wc->wcc[i].n, cyc_sum[i], tot);
742 /* Print timing information when it is for a PP or PP+PME
744 print_cycles(fplog, c2t_pp, wcn[i],
746 wc->wcc[i].n, cyc_sum[i], tot);
747 tot_for_pp += cyc_sum[i];
750 if (wc->wcc_all != NULL)
752 for (i = 0; i < ewcNR; i++)
754 for (j = 0; j < ewcNR; j++)
756 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
757 print_cycles(fplog, c2t_pp, buf,
759 wc->wcc_all[i*ewcNR+j].n,
760 wc->wcc_all[i*ewcNR+j].c,
765 tot_for_rest = tot * npp * nth_pp / static_cast<double>(nth_tot);
766 print_cycles(fplog, c2t_pp, "Rest",
768 -1, tot_for_rest - tot_for_pp, tot);
769 fprintf(fplog, "%s\n", hline);
770 print_cycles(fplog, c2t, "Total",
773 fprintf(fplog, "%s\n", hline);
778 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
779 " twice the total reported, but the cycle count total and %% are correct.\n"
783 if (wc->wcc[ewcPMEMESH].n > 0)
785 fprintf(fplog, " Breakdown of PME mesh computation\n");
786 fprintf(fplog, "%s\n", hline);
787 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
789 if (is_pme_subcounter(i))
791 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i],
792 npme > 0 ? npme : npp, nth_pme,
793 wc->wcc[i].n, cyc_sum[i], tot);
796 fprintf(fplog, "%s\n", hline);
799 #ifdef GMX_CYCLE_SUBCOUNTERS
800 fprintf(fplog, " Breakdown of PP computation\n");
801 fprintf(fplog, "%s\n", hline);
802 for (i = 0; i < ewcsNR; i++)
804 print_cycles(fplog, c2t_pp, wcsn[i],
806 wc->wcsc[i].n, cyc_sum[ewcNR+i], tot);
808 fprintf(fplog, "%s\n", hline);
811 /* print GPU timing summary */
814 const char *k_log_str[2][2] = {
815 {"Nonbonded F kernel", "Nonbonded F+ene k."},
816 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
819 tot_gpu = gpu_t->pl_h2d_t + gpu_t->nb_h2d_t + gpu_t->nb_d2h_t;
821 /* add up the kernel timings */
823 for (i = 0; i < 2; i++)
825 for (j = 0; j < 2; j++)
827 tot_k += gpu_t->ktime[i][j].t;
832 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
833 if (wc->wcc[ewcPMEMESH].n > 0)
835 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
837 tot_cpu_overlap *= realtime*1000/tot; /* convert s to ms */
839 fprintf(fplog, "\n GPU timings\n%s\n", hline);
840 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
841 fprintf(fplog, "%s\n", hline);
842 print_gputimes(fplog, "Pair list H2D",
843 gpu_t->pl_h2d_c, gpu_t->pl_h2d_t, tot_gpu);
844 print_gputimes(fplog, "X / q H2D",
845 gpu_t->nb_c, gpu_t->nb_h2d_t, tot_gpu);
847 for (i = 0; i < 2; i++)
849 for (j = 0; j < 2; j++)
851 if (gpu_t->ktime[i][j].c)
853 print_gputimes(fplog, k_log_str[i][j],
854 gpu_t->ktime[i][j].c, gpu_t->ktime[i][j].t, tot_gpu);
859 print_gputimes(fplog, "F D2H", gpu_t->nb_c, gpu_t->nb_d2h_t, tot_gpu);
860 fprintf(fplog, "%s\n", hline);
861 print_gputimes(fplog, "Total ", gpu_t->nb_c, tot_gpu, tot_gpu);
862 fprintf(fplog, "%s\n", hline);
864 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
865 if (gpu_t->nb_c > 0 && wc->wcc[ewcFORCE].n > 0)
867 fprintf(fplog, "\nForce evaluation time GPU/CPU: %.3f ms/%.3f ms = %.3f\n",
868 tot_gpu/gpu_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
872 /* only print notes related to CPU-GPU load balance with PME */
873 if (wc->wcc[ewcPMEMESH].n > 0)
875 fprintf(fplog, "For optimal performance this ratio should be close to 1!\n");
877 /* print note if the imbalance is high with PME case in which
878 * CPU-GPU load balancing is possible */
879 if (gpu_cpu_ratio < 0.75 || gpu_cpu_ratio > 1.2)
881 /* Only the sim master calls this function, so always print to stderr */
882 if (gpu_cpu_ratio < 0.75)
886 /* The user could have used -notunepme,
887 * but we currently can't check that here.
889 md_print_warn(NULL, fplog,
890 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
891 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
892 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
896 /* We should not end up here, unless the box is
897 * too small for increasing the cut-off for PME tuning.
899 md_print_warn(NULL, fplog,
900 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
901 " performance loss.");
904 if (gpu_cpu_ratio > 1.2)
906 md_print_warn(NULL, fplog,
907 "\nNOTE: The GPU has >20%% more load than the CPU. This imbalance causes\n"
908 " performance loss, consider using a shorter cut-off and a finer PME grid.");
916 md_print_warn(NULL, fplog,
917 "MPI_Barrier was called before each cycle start/stop\n"
918 "call, so timings are not those of real runs.\n");
921 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
922 (cyc_sum[ewcDOMDEC] > tot*0.1 ||
923 cyc_sum[ewcNS] > tot*0.1))
925 /* Only the sim master calls this function, so always print to stderr */
926 if (wc->wcc[ewcDOMDEC].n == 0)
928 md_print_warn(NULL, fplog,
929 "NOTE: %d %% of the run time was spent in pair search,\n"
930 " you might want to increase nstlist (this has no effect on accuracy)\n",
931 (int)(100*cyc_sum[ewcNS]/tot+0.5));
935 md_print_warn(NULL, fplog,
936 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
937 " %d %% of the run time was spent in pair search,\n"
938 " you might want to increase nstlist (this has no effect on accuracy)\n",
939 (int)(100*cyc_sum[ewcDOMDEC]/tot+0.5),
940 (int)(100*cyc_sum[ewcNS]/tot+0.5));
944 if (cyc_sum[ewcMoveE] > tot*0.05)
946 /* Only the sim master calls this function, so always print to stderr */
947 md_print_warn(NULL, fplog,
948 "NOTE: %d %% of the run time was spent communicating energies,\n"
949 " you might want to use the -gcom option of mdrun\n",
950 (int)(100*cyc_sum[ewcMoveE]/tot+0.5));
954 extern gmx_int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
961 return wc->reset_counters;
964 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_int64_t reset_counters)
971 wc->reset_counters = reset_counters;
974 #ifdef GMX_CYCLE_SUBCOUNTERS
976 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
980 wc->wcsc[ewcs].start = gmx_cycles_read();
984 void wallcycle_sub_start_nocount(gmx_wallcycle_t wc, int ewcs)
991 wallcycle_sub_start(wc, ewcs);
995 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
999 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
1006 void wallcycle_sub_start(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
1009 void wallcycle_sub_start_nocount(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
1012 void wallcycle_sub_stop(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
1016 #endif /* GMX_CYCLE_SUBCOUNTERS */