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, 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/utility/cstringutil.h"
49 #include "gromacs/utility/gmxmpi.h"
50 #include "gromacs/utility/smalloc.h"
52 /* DEBUG_WCYCLE adds consistency checking for the counters.
53 * It checks if you stop a counter different from the last
54 * one that was opened and if you do nest too deep.
56 /* #define DEBUG_WCYCLE */
59 #include "gromacs/utility/fatalerror.h"
70 typedef struct gmx_wallcycle
73 /* variables for testing/debugging */
79 int counterlist[DEPTH_MAX];
83 gmx_cycles_t cycle_prev;
84 gmx_int64_t reset_counters;
86 MPI_Comm mpi_comm_mygroup;
90 #ifdef GMX_CYCLE_SUBCOUNTERS
96 /* Each name should not exceed 19 printing characters
97 (ie. terminating null can be twentieth) */
98 static const char *wcn[ewcNR] =
100 "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
101 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
102 "Comm. coord.", "Born radii", "Force", "Wait + Comm. F", "PME mesh",
103 "PME redist. X/F", "PME spread/gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve LJ", "PME solve Elec",
104 "PME wait for PP", "Wait + Recv. PME F", "Wait GPU nonlocal", "Wait GPU local", "Wait GPU loc. est.", "NB X/F buffer ops.",
105 "Vsite spread", "COM pull force",
106 "Write traj.", "Update", "Constraints", "Comm. energies",
107 "Enforced rotation", "Add rot. forces", "Coordinate swapping", "IMD", "Test"
110 static const char *wcsn[ewcsNR] =
112 "DD redist.", "DD NS grid + sort", "DD setup comm.",
113 "DD make top.", "DD make constr.", "DD top. other",
114 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
115 "Listed F", "Nonbonded F", "Ewald F correction",
116 "NB X buffer ops.", "NB F buffer ops."
119 gmx_bool wallcycle_have_counter(void)
121 return gmx_cycles_have_counter();
124 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec gmx_unused *cr,
125 int nthreads_pp, int nthreads_pme)
130 if (!wallcycle_have_counter())
137 wc->wc_barrier = FALSE;
141 wc->reset_counters = resetstep;
142 wc->nthreads_pp = nthreads_pp;
143 wc->nthreads_pme = nthreads_pme;
144 wc->cycles_sum = NULL;
147 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != NULL)
151 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
153 wc->wc_barrier = TRUE;
154 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
158 snew(wc->wcc, ewcNR);
159 if (getenv("GMX_CYCLE_ALL") != NULL)
163 fprintf(fplog, "\nWill time all the code during the run\n\n");
165 snew(wc->wcc_all, ewcNR*ewcNR);
168 #ifdef GMX_CYCLE_SUBCOUNTERS
169 snew(wc->wcsc, ewcsNR);
179 void wallcycle_destroy(gmx_wallcycle_t wc)
190 if (wc->wcc_all != NULL)
194 #ifdef GMX_CYCLE_SUBCOUNTERS
195 if (wc->wcsc != NULL)
203 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
206 wc->cycle_prev = cycle;
209 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
211 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
212 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
217 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
219 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
221 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
223 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
226 wc->counterlist[wc->count_depth] = ewc;
230 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
234 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
236 if (wc->count_depth < 0)
238 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
240 if (wc->counterlist[wc->count_depth] != ewc)
242 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
243 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
248 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
260 MPI_Barrier(wc->mpi_comm_mygroup);
265 debug_start_check(wc, ewc);
268 cycle = gmx_cycles_read();
269 wc->wcc[ewc].start = cycle;
270 if (wc->wcc_all != NULL)
275 wallcycle_all_start(wc, ewc, cycle);
277 else if (wc->wc_depth == 3)
279 wallcycle_all_stop(wc, ewc, cycle);
284 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
291 wallcycle_start(wc, ewc);
295 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
297 gmx_cycles_t cycle, last;
307 MPI_Barrier(wc->mpi_comm_mygroup);
312 debug_stop_check(wc, ewc);
315 cycle = gmx_cycles_read();
316 last = cycle - wc->wcc[ewc].start;
317 wc->wcc[ewc].c += last;
324 wallcycle_all_stop(wc, ewc, cycle);
326 else if (wc->wc_depth == 2)
328 wallcycle_all_start(wc, ewc, cycle);
335 void wallcycle_reset_all(gmx_wallcycle_t wc)
344 for (i = 0; i < ewcNR; i++)
351 for (i = 0; i < ewcNR*ewcNR; i++)
353 wc->wcc_all[i].n = 0;
354 wc->wcc_all[i].c = 0;
357 #ifdef GMX_CYCLE_SUBCOUNTERS
358 for (i = 0; i < ewcsNR; i++)
366 static gmx_bool is_pme_counter(int ewc)
368 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
371 static gmx_bool is_pme_subcounter(int ewc)
373 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
376 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
379 double cycles[ewcNR+ewcsNR];
380 double cycles_n[ewcNR+ewcsNR], buf[ewcNR+ewcsNR], *cyc_all, *buf_all;
389 snew(wc->cycles_sum, ewcNR+ewcsNR);
393 /* The GPU wait estimate counter is used for load balancing only
394 * and will mess up the total due to double counting: clear it.
396 wcc[ewcWAIT_GPU_NB_L_EST].n = 0;
397 wcc[ewcWAIT_GPU_NB_L_EST].c = 0;
399 for (i = 0; i < ewcNR; i++)
401 if (is_pme_counter(i) || (i == ewcRUN && cr->duty == DUTY_PME))
403 wcc[i].c *= wc->nthreads_pme;
407 for (j = 0; j < ewcNR; j++)
409 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pme;
415 wcc[i].c *= wc->nthreads_pp;
419 for (j = 0; j < ewcNR; j++)
421 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pp;
427 if (wcc[ewcDDCOMMLOAD].n > 0)
429 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMLOAD].c;
431 if (wcc[ewcDDCOMMBOUND].n > 0)
433 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMBOUND].c;
435 if (wcc[ewcPME_FFTCOMM].n > 0)
437 wcc[ewcPME_FFT].c -= wcc[ewcPME_FFTCOMM].c;
440 if (cr->npmenodes == 0)
442 /* All nodes do PME (or no PME at all) */
443 if (wcc[ewcPMEMESH].n > 0)
445 wcc[ewcFORCE].c -= wcc[ewcPMEMESH].c;
450 /* The are PME-only nodes */
451 if (wcc[ewcPMEMESH].n > 0)
453 /* This must be a PME only node, calculate the Wait + Comm. time */
454 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
458 /* Store the cycles in a double buffer for summing */
459 for (i = 0; i < ewcNR; i++)
461 cycles_n[i] = (double)wcc[i].n;
462 cycles[i] = (double)wcc[i].c;
465 #ifdef GMX_CYCLE_SUBCOUNTERS
466 for (i = 0; i < ewcsNR; i++)
468 wc->wcsc[i].c *= wc->nthreads_pp;
469 cycles_n[ewcNR+i] = (double)wc->wcsc[i].n;
470 cycles[ewcNR+i] = (double)wc->wcsc[i].c;
478 MPI_Allreduce(cycles_n, buf, nsum, MPI_DOUBLE, MPI_MAX,
480 for (i = 0; i < ewcNR; i++)
482 wcc[i].n = (int)(buf[i] + 0.5);
484 #ifdef GMX_CYCLE_SUBCOUNTERS
485 for (i = 0; i < ewcsNR; i++)
487 wc->wcsc[i].n = (int)(buf[ewcNR+i] + 0.5);
491 MPI_Allreduce(cycles, wc->cycles_sum, nsum, MPI_DOUBLE, MPI_SUM,
494 if (wc->wcc_all != NULL)
496 snew(cyc_all, ewcNR*ewcNR);
497 snew(buf_all, ewcNR*ewcNR);
498 for (i = 0; i < ewcNR*ewcNR; i++)
500 cyc_all[i] = wc->wcc_all[i].c;
502 MPI_Allreduce(cyc_all, buf_all, ewcNR*ewcNR, MPI_DOUBLE, MPI_SUM,
504 for (i = 0; i < ewcNR*ewcNR; i++)
506 wc->wcc_all[i].c = buf_all[i];
515 for (i = 0; i < nsum; i++)
517 wc->cycles_sum[i] = cycles[i];
522 static void print_cycles(FILE *fplog, double c2t, const char *name,
523 int nnodes, int nthreads,
524 int ncalls, double c_sum, double tot)
527 char nthreads_str[6];
535 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
538 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
542 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
546 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
550 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
559 /* Convert the cycle count to wallclock time for this task */
562 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
563 name, nnodes_str, nthreads_str, ncalls_str, wallt,
564 c_sum*1e-9, 100*c_sum/tot);
568 static void print_gputimes(FILE *fplog, const char *name,
569 int n, double t, double tot_t)
576 snprintf(num, sizeof(num), "%10d", n);
577 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
582 sprintf(avg_perf, " ");
586 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
587 name, num, t/1000, avg_perf, 100 * t/tot_t);
591 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
592 name, "", t/1000, avg_perf, 100.0);
596 static void print_header(FILE *fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
598 int nrank_tot = nrank_pp + nrank_pme;
601 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
604 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
606 /* Don't report doing PP+PME, because we can't tell here if
607 * this is RF, etc. */
611 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
614 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
616 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
619 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
623 fprintf(fplog, "\n\n");
624 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
625 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
628 void wallcycle_print(FILE *fplog, int nnodes, int npme, double realtime,
629 gmx_wallcycle_t wc, wallclock_gpu_t *gpu_t)
632 double tot, tot_for_pp, tot_for_rest, tot_gpu, tot_cpu_overlap, gpu_cpu_ratio, tot_k;
633 double c2t, c2t_pp, c2t_pme = 0;
634 int i, j, npp, nth_pp, nth_pme, nth_tot;
636 const char *hline = "-----------------------------------------------------------------------------";
643 nth_pp = wc->nthreads_pp;
644 nth_pme = wc->nthreads_pme;
646 cyc_sum = wc->cycles_sum;
649 nth_tot = npp*nth_pp + npme*nth_pme;
651 /* When using PME-only nodes, the next line is valid for both
652 PP-only and PME-only nodes because they started ewcRUN at the
654 tot = cyc_sum[ewcRUN];
657 /* Conversion factor from cycles to seconds */
661 c2t_pp = c2t * nth_tot / (double) (npp*nth_pp);
664 c2t_pme = c2t * nth_tot / (double) (npme*nth_pme);
674 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");
676 print_header(fplog, npp, nth_pp, npme, nth_pme);
678 fprintf(fplog, "%s\n", hline);
679 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
681 if (is_pme_subcounter(i))
683 /* Do not count these at all */
685 else if (npme > 0 && is_pme_counter(i))
687 /* Print timing information for PME-only nodes, but add an
688 * asterisk so the reader of the table can know that the
689 * walltimes are not meant to add up. The asterisk still
690 * fits in the required maximum of 19 characters. */
692 snprintf(buffer, STRLEN, "%s *", wcn[i]);
693 print_cycles(fplog, c2t_pme, buffer,
695 wc->wcc[i].n, cyc_sum[i], tot);
699 /* Print timing information when it is for a PP or PP+PME
701 print_cycles(fplog, c2t_pp, wcn[i],
703 wc->wcc[i].n, cyc_sum[i], tot);
704 tot_for_pp += cyc_sum[i];
707 if (wc->wcc_all != NULL)
709 for (i = 0; i < ewcNR; i++)
711 for (j = 0; j < ewcNR; j++)
713 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
714 print_cycles(fplog, c2t_pp, buf,
716 wc->wcc_all[i*ewcNR+j].n,
717 wc->wcc_all[i*ewcNR+j].c,
722 tot_for_rest = tot * (npp * nth_pp) / (double) nth_tot;
723 print_cycles(fplog, c2t_pp, "Rest",
725 -1, tot_for_rest - tot_for_pp, tot);
726 fprintf(fplog, "%s\n", hline);
727 print_cycles(fplog, c2t, "Total",
730 fprintf(fplog, "%s\n", hline);
735 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
736 " twice the total reported, but the cycle count total and %% are correct.\n"
740 if (wc->wcc[ewcPMEMESH].n > 0)
742 fprintf(fplog, " Breakdown of PME mesh computation\n");
743 fprintf(fplog, "%s\n", hline);
744 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
746 if (is_pme_subcounter(i))
748 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i],
749 npme > 0 ? npme : npp, nth_pme,
750 wc->wcc[i].n, cyc_sum[i], tot);
753 fprintf(fplog, "%s\n", hline);
756 #ifdef GMX_CYCLE_SUBCOUNTERS
757 fprintf(fplog, " Breakdown of PP computation\n");
758 fprintf(fplog, "%s\n", hline);
759 for (i = 0; i < ewcsNR; i++)
761 print_cycles(fplog, c2t_pp, wcsn[i],
763 wc->wcsc[i].n, cyc_sum[ewcNR+i], tot);
765 fprintf(fplog, "%s\n", hline);
768 /* print GPU timing summary */
771 const char *k_log_str[2][2] = {
772 {"Nonbonded F kernel", "Nonbonded F+ene k."},
773 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
776 tot_gpu = gpu_t->pl_h2d_t + gpu_t->nb_h2d_t + gpu_t->nb_d2h_t;
778 /* add up the kernel timings */
780 for (i = 0; i < 2; i++)
782 for (j = 0; j < 2; j++)
784 tot_k += gpu_t->ktime[i][j].t;
789 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
790 if (wc->wcc[ewcPMEMESH].n > 0)
792 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
794 tot_cpu_overlap *= realtime*1000/tot; /* convert s to ms */
796 fprintf(fplog, "\n GPU timings\n%s\n", hline);
797 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
798 fprintf(fplog, "%s\n", hline);
799 print_gputimes(fplog, "Pair list H2D",
800 gpu_t->pl_h2d_c, gpu_t->pl_h2d_t, tot_gpu);
801 print_gputimes(fplog, "X / q H2D",
802 gpu_t->nb_c, gpu_t->nb_h2d_t, tot_gpu);
804 for (i = 0; i < 2; i++)
806 for (j = 0; j < 2; j++)
808 if (gpu_t->ktime[i][j].c)
810 print_gputimes(fplog, k_log_str[i][j],
811 gpu_t->ktime[i][j].c, gpu_t->ktime[i][j].t, tot_gpu);
816 print_gputimes(fplog, "F D2H", gpu_t->nb_c, gpu_t->nb_d2h_t, tot_gpu);
817 fprintf(fplog, "%s\n", hline);
818 print_gputimes(fplog, "Total ", gpu_t->nb_c, tot_gpu, tot_gpu);
819 fprintf(fplog, "%s\n", hline);
821 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
822 fprintf(fplog, "\nForce evaluation time GPU/CPU: %.3f ms/%.3f ms = %.3f\n",
823 tot_gpu/gpu_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
826 /* only print notes related to CPU-GPU load balance with PME */
827 if (wc->wcc[ewcPMEMESH].n > 0)
829 fprintf(fplog, "For optimal performance this ratio should be close to 1!\n");
831 /* print note if the imbalance is high with PME case in which
832 * CPU-GPU load balancing is possible */
833 if (gpu_cpu_ratio < 0.75 || gpu_cpu_ratio > 1.2)
835 /* Only the sim master calls this function, so always print to stderr */
836 if (gpu_cpu_ratio < 0.75)
840 /* The user could have used -notunepme,
841 * but we currently can't check that here.
843 md_print_warn(NULL, fplog,
844 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
845 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
846 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
850 /* We should not end up here, unless the box is
851 * too small for increasing the cut-off for PME tuning.
853 md_print_warn(NULL, fplog,
854 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
855 " performance loss.");
858 if (gpu_cpu_ratio > 1.2)
860 md_print_warn(NULL, fplog,
861 "\nNOTE: The GPU has >20%% more load than the CPU. This imbalance causes\n"
862 " performance loss, consider using a shorter cut-off and a finer PME grid.");
868 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
869 (cyc_sum[ewcDOMDEC] > tot*0.1 ||
870 cyc_sum[ewcNS] > tot*0.1))
872 /* Only the sim master calls this function, so always print to stderr */
873 if (wc->wcc[ewcDOMDEC].n == 0)
875 md_print_warn(NULL, fplog,
876 "NOTE: %d %% of the run time was spent in pair search,\n"
877 " you might want to increase nstlist (this has no effect on accuracy)\n",
878 (int)(100*cyc_sum[ewcNS]/tot+0.5));
882 md_print_warn(NULL, fplog,
883 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
884 " %d %% of the run time was spent in pair search,\n"
885 " you might want to increase nstlist (this has no effect on accuracy)\n",
886 (int)(100*cyc_sum[ewcDOMDEC]/tot+0.5),
887 (int)(100*cyc_sum[ewcNS]/tot+0.5));
891 if (cyc_sum[ewcMoveE] > tot*0.05)
893 /* Only the sim master calls this function, so always print to stderr */
894 md_print_warn(NULL, fplog,
895 "NOTE: %d %% of the run time was spent communicating energies,\n"
896 " you might want to use the -gcom option of mdrun\n",
897 (int)(100*cyc_sum[ewcMoveE]/tot+0.5));
901 extern gmx_int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
908 return wc->reset_counters;
911 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_int64_t reset_counters)
918 wc->reset_counters = reset_counters;
921 #ifdef GMX_CYCLE_SUBCOUNTERS
923 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
927 wc->wcsc[ewcs].start = gmx_cycles_read();
931 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
935 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
942 void wallcycle_sub_start(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
945 void wallcycle_sub_stop(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
949 #endif /* GMX_CYCLE_SUBCOUNTERS */