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.
37 #include "gromacs/timing/wallcycle.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "gmx_fatal.h"
47 #include "md_logging.h"
48 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/timing/cyclecounter.h"
51 #include "gromacs/utility/gmxmpi.h"
53 /* DEBUG_WCYCLE adds consistency checking for the counters.
54 * It checks if you stop a counter different from the last
55 * one that was opened and if you do nest too deep.
57 /* #define DEBUG_WCYCLE */
67 typedef struct gmx_wallcycle
70 /* variables for testing/debugging */
76 int counterlist[DEPTH_MAX];
80 gmx_cycles_t cycle_prev;
81 gmx_int64_t reset_counters;
83 MPI_Comm mpi_comm_mygroup;
87 #ifdef GMX_CYCLE_SUBCOUNTERS
93 /* Each name should not exceed 19 printing characters
94 (ie. terminating null can be twentieth) */
95 static const char *wcn[ewcNR] =
97 "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
98 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
99 "Comm. coord.", "Born radii", "Force", "Wait + Comm. F", "PME mesh",
100 "PME redist. X/F", "PME spread/gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve LJ", "PME solve Elec",
101 "PME wait for PP", "Wait + Recv. PME F", "Wait GPU nonlocal", "Wait GPU local", "NB X/F buffer ops.",
102 "Vsite spread", "Write traj.", "Update", "Constraints", "Comm. energies",
103 "Enforced rotation", "Add rot. forces", "Coordinate swapping", "IMD", "Test"
106 static const char *wcsn[ewcsNR] =
108 "DD redist.", "DD NS grid + sort", "DD setup comm.",
109 "DD make top.", "DD make constr.", "DD top. other",
110 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
111 "Bonded F", "Nonbonded F", "Ewald F correction",
112 "NB X buffer ops.", "NB F buffer ops."
115 gmx_bool wallcycle_have_counter(void)
117 return gmx_cycles_have_counter();
120 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec gmx_unused *cr,
121 int nthreads_pp, int nthreads_pme)
126 if (!wallcycle_have_counter())
133 wc->wc_barrier = FALSE;
137 wc->reset_counters = resetstep;
138 wc->nthreads_pp = nthreads_pp;
139 wc->nthreads_pme = nthreads_pme;
140 wc->cycles_sum = NULL;
143 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != NULL)
147 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
149 wc->wc_barrier = TRUE;
150 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
154 snew(wc->wcc, ewcNR);
155 if (getenv("GMX_CYCLE_ALL") != NULL)
159 fprintf(fplog, "\nWill time all the code during the run\n\n");
161 snew(wc->wcc_all, ewcNR*ewcNR);
164 #ifdef GMX_CYCLE_SUBCOUNTERS
165 snew(wc->wcsc, ewcsNR);
175 void wallcycle_destroy(gmx_wallcycle_t wc)
186 if (wc->wcc_all != NULL)
190 #ifdef GMX_CYCLE_SUBCOUNTERS
191 if (wc->wcsc != NULL)
199 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
202 wc->cycle_prev = cycle;
205 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
207 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
208 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
213 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
215 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
217 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
219 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
222 wc->counterlist[wc->count_depth] = ewc;
226 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
230 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
232 if (wc->count_depth < 0)
234 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
236 if (wc->counterlist[wc->count_depth] != ewc)
238 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
239 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
244 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
256 MPI_Barrier(wc->mpi_comm_mygroup);
261 debug_start_check(wc, ewc);
264 cycle = gmx_cycles_read();
265 wc->wcc[ewc].start = cycle;
266 if (wc->wcc_all != NULL)
271 wallcycle_all_start(wc, ewc, cycle);
273 else if (wc->wc_depth == 3)
275 wallcycle_all_stop(wc, ewc, cycle);
280 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
287 wallcycle_start(wc, ewc);
291 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
293 gmx_cycles_t cycle, last;
303 MPI_Barrier(wc->mpi_comm_mygroup);
308 debug_stop_check(wc, ewc);
311 cycle = gmx_cycles_read();
312 last = cycle - wc->wcc[ewc].start;
313 wc->wcc[ewc].c += last;
320 wallcycle_all_stop(wc, ewc, cycle);
322 else if (wc->wc_depth == 2)
324 wallcycle_all_start(wc, ewc, cycle);
331 void wallcycle_reset_all(gmx_wallcycle_t wc)
340 for (i = 0; i < ewcNR; i++)
347 for (i = 0; i < ewcNR*ewcNR; i++)
349 wc->wcc_all[i].n = 0;
350 wc->wcc_all[i].c = 0;
353 #ifdef GMX_CYCLE_SUBCOUNTERS
354 for (i = 0; i < ewcsNR; i++)
362 static gmx_bool is_pme_counter(int ewc)
364 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
367 static gmx_bool is_pme_subcounter(int ewc)
369 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
372 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
375 double cycles[ewcNR+ewcsNR];
376 double cycles_n[ewcNR+ewcsNR], buf[ewcNR+ewcsNR], *cyc_all, *buf_all;
385 snew(wc->cycles_sum, ewcNR+ewcsNR);
389 for (i = 0; i < ewcNR; i++)
391 if (is_pme_counter(i) || (i == ewcRUN && cr->duty == DUTY_PME))
393 wcc[i].c *= wc->nthreads_pme;
397 for (j = 0; j < ewcNR; j++)
399 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pme;
405 wcc[i].c *= wc->nthreads_pp;
409 for (j = 0; j < ewcNR; j++)
411 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pp;
417 if (wcc[ewcDDCOMMLOAD].n > 0)
419 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMLOAD].c;
421 if (wcc[ewcDDCOMMBOUND].n > 0)
423 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMBOUND].c;
425 if (wcc[ewcPME_FFTCOMM].n > 0)
427 wcc[ewcPME_FFT].c -= wcc[ewcPME_FFTCOMM].c;
430 if (cr->npmenodes == 0)
432 /* All nodes do PME (or no PME at all) */
433 if (wcc[ewcPMEMESH].n > 0)
435 wcc[ewcFORCE].c -= wcc[ewcPMEMESH].c;
440 /* The are PME-only nodes */
441 if (wcc[ewcPMEMESH].n > 0)
443 /* This must be a PME only node, calculate the Wait + Comm. time */
444 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
448 /* Store the cycles in a double buffer for summing */
449 for (i = 0; i < ewcNR; i++)
451 cycles_n[i] = (double)wcc[i].n;
452 cycles[i] = (double)wcc[i].c;
455 #ifdef GMX_CYCLE_SUBCOUNTERS
456 for (i = 0; i < ewcsNR; i++)
458 wc->wcsc[i].c *= wc->nthreads_pp;
459 cycles_n[ewcNR+i] = (double)wc->wcsc[i].n;
460 cycles[ewcNR+i] = (double)wc->wcsc[i].c;
468 MPI_Allreduce(cycles_n, buf, nsum, MPI_DOUBLE, MPI_MAX,
470 for (i = 0; i < ewcNR; i++)
472 wcc[i].n = (int)(buf[i] + 0.5);
474 #ifdef GMX_CYCLE_SUBCOUNTERS
475 for (i = 0; i < ewcsNR; i++)
477 wc->wcsc[i].n = (int)(buf[ewcNR+i] + 0.5);
481 MPI_Allreduce(cycles, wc->cycles_sum, nsum, MPI_DOUBLE, MPI_SUM,
484 if (wc->wcc_all != NULL)
486 snew(cyc_all, ewcNR*ewcNR);
487 snew(buf_all, ewcNR*ewcNR);
488 for (i = 0; i < ewcNR*ewcNR; i++)
490 cyc_all[i] = wc->wcc_all[i].c;
492 MPI_Allreduce(cyc_all, buf_all, ewcNR*ewcNR, MPI_DOUBLE, MPI_SUM,
494 for (i = 0; i < ewcNR*ewcNR; i++)
496 wc->wcc_all[i].c = buf_all[i];
505 for (i = 0; i < nsum; i++)
507 wc->cycles_sum[i] = cycles[i];
512 static void print_cycles(FILE *fplog, double c2t, const char *name,
513 int nnodes, int nthreads,
514 int ncalls, double c_sum, double tot)
517 char nthreads_str[6];
525 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
528 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
532 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
536 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
540 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
549 /* Convert the cycle count to wallclock time for this task */
552 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
553 name, nnodes_str, nthreads_str, ncalls_str, wallt,
554 c_sum*1e-9, 100*c_sum/tot);
558 static void print_gputimes(FILE *fplog, const char *name,
559 int n, double t, double tot_t)
566 snprintf(num, sizeof(num), "%10d", n);
567 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
572 sprintf(avg_perf, " ");
576 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
577 name, num, t/1000, avg_perf, 100 * t/tot_t);
581 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
582 name, "", t/1000, avg_perf, 100.0);
586 static void print_header(FILE *fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
588 int nrank_tot = nrank_pp + nrank_pme;
591 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
594 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
596 /* Don't report doing PP+PME, because we can't tell here if
597 * this is RF, etc. */
601 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
604 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
606 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
609 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
613 fprintf(fplog, "\n\n");
614 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
615 fprintf(fplog, " Nodes Threads Count (s) total sum %%\n");
618 void wallcycle_print(FILE *fplog, int nnodes, int npme, double realtime,
619 gmx_wallcycle_t wc, wallclock_gpu_t *gpu_t)
622 double tot, tot_for_pp, tot_for_rest, tot_gpu, tot_cpu_overlap, gpu_cpu_ratio, tot_k;
623 double c2t, c2t_pp, c2t_pme;
624 int i, j, npp, nth_pp, nth_pme, nth_tot;
626 const char *hline = "-----------------------------------------------------------------------------";
633 nth_pp = wc->nthreads_pp;
634 nth_pme = wc->nthreads_pme;
636 cyc_sum = wc->cycles_sum;
639 nth_tot = npp*nth_pp + npme*nth_pme;
641 /* When using PME-only nodes, the next line is valid for both
642 PP-only and PME-only nodes because they started ewcRUN at the
644 tot = cyc_sum[ewcRUN];
647 /* Conversion factor from cycles to seconds */
651 c2t_pp = c2t * nth_tot / (double) (npp*nth_pp);
652 c2t_pme = c2t * nth_tot / (double) (npme*nth_pme);
661 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");
663 print_header(fplog, npp, nth_pp, npme, nth_pme);
665 fprintf(fplog, "%s\n", hline);
666 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
668 if (is_pme_subcounter(i))
670 /* Do not count these at all */
672 else if (npme > 0 && is_pme_counter(i))
674 /* Print timing information for PME-only nodes, but add an
675 * asterisk so the reader of the table can know that the
676 * walltimes are not meant to add up. The asterisk still
677 * fits in the required maximum of 19 characters. */
679 snprintf(buffer, STRLEN, "%s *", wcn[i]);
680 print_cycles(fplog, c2t_pme, buffer,
682 wc->wcc[i].n, cyc_sum[i], tot);
686 /* Print timing information when it is for a PP or PP+PME
688 print_cycles(fplog, c2t_pp, wcn[i],
690 wc->wcc[i].n, cyc_sum[i], tot);
691 tot_for_pp += cyc_sum[i];
694 if (wc->wcc_all != NULL)
696 for (i = 0; i < ewcNR; i++)
698 for (j = 0; j < ewcNR; j++)
700 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
701 print_cycles(fplog, c2t_pp, buf,
703 wc->wcc_all[i*ewcNR+j].n,
704 wc->wcc_all[i*ewcNR+j].c,
709 tot_for_rest = tot * (npp * nth_pp) / (double) nth_tot;
710 print_cycles(fplog, c2t_pp, "Rest",
712 -1, tot_for_rest - tot_for_pp, tot);
713 fprintf(fplog, "%s\n", hline);
714 print_cycles(fplog, c2t, "Total",
717 fprintf(fplog, "%s\n", hline);
722 "(*) Note that with separate PME nodes, the walltime column actually sums to\n"
723 " twice the total reported, but the cycle count total and %% are correct.\n"
727 if (wc->wcc[ewcPMEMESH].n > 0)
729 fprintf(fplog, " Breakdown of PME mesh computation\n");
730 fprintf(fplog, "%s\n", hline);
731 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
733 if (is_pme_subcounter(i))
735 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i],
736 npme > 0 ? npme : npp, nth_pme,
737 wc->wcc[i].n, cyc_sum[i], tot);
740 fprintf(fplog, "%s\n", hline);
743 #ifdef GMX_CYCLE_SUBCOUNTERS
744 fprintf(fplog, " Breakdown of PP computation\n");
745 fprintf(fplog, "%s\n", hline);
746 for (i = 0; i < ewcsNR; i++)
748 print_cycles(fplog, c2t_pp, wcsn[i],
750 wc->wcsc[i].n, cyc_sum[ewcNR+i], tot);
752 fprintf(fplog, "%s\n", hline);
755 /* print GPU timing summary */
758 const char *k_log_str[2][2] = {
759 {"Nonbonded F kernel", "Nonbonded F+ene k."},
760 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
763 tot_gpu = gpu_t->pl_h2d_t + gpu_t->nb_h2d_t + gpu_t->nb_d2h_t;
765 /* add up the kernel timings */
767 for (i = 0; i < 2; i++)
769 for (j = 0; j < 2; j++)
771 tot_k += gpu_t->ktime[i][j].t;
776 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
777 if (wc->wcc[ewcPMEMESH].n > 0)
779 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
781 tot_cpu_overlap *= realtime*1000/tot; /* convert s to ms */
783 fprintf(fplog, "\n GPU timings\n%s\n", hline);
784 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
785 fprintf(fplog, "%s\n", hline);
786 print_gputimes(fplog, "Pair list H2D",
787 gpu_t->pl_h2d_c, gpu_t->pl_h2d_t, tot_gpu);
788 print_gputimes(fplog, "X / q H2D",
789 gpu_t->nb_c, gpu_t->nb_h2d_t, tot_gpu);
791 for (i = 0; i < 2; i++)
793 for (j = 0; j < 2; j++)
795 if (gpu_t->ktime[i][j].c)
797 print_gputimes(fplog, k_log_str[i][j],
798 gpu_t->ktime[i][j].c, gpu_t->ktime[i][j].t, tot_gpu);
803 print_gputimes(fplog, "F D2H", gpu_t->nb_c, gpu_t->nb_d2h_t, tot_gpu);
804 fprintf(fplog, "%s\n", hline);
805 print_gputimes(fplog, "Total ", gpu_t->nb_c, tot_gpu, tot_gpu);
806 fprintf(fplog, "%s\n", hline);
808 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
809 fprintf(fplog, "\nForce evaluation time GPU/CPU: %.3f ms/%.3f ms = %.3f\n",
810 tot_gpu/gpu_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
813 /* only print notes related to CPU-GPU load balance with PME */
814 if (wc->wcc[ewcPMEMESH].n > 0)
816 fprintf(fplog, "For optimal performance this ratio should be close to 1!\n");
818 /* print note if the imbalance is high with PME case in which
819 * CPU-GPU load balancing is possible */
820 if (gpu_cpu_ratio < 0.75 || gpu_cpu_ratio > 1.2)
822 /* Only the sim master calls this function, so always print to stderr */
823 if (gpu_cpu_ratio < 0.75)
827 /* The user could have used -notunepme,
828 * but we currently can't check that here.
830 md_print_warn(NULL, fplog,
831 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
832 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
833 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
837 /* We should not end up here, unless the box is
838 * too small for increasing the cut-off for PME tuning.
840 md_print_warn(NULL, fplog,
841 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
842 " performance loss.");
845 if (gpu_cpu_ratio > 1.2)
847 md_print_warn(NULL, fplog,
848 "\nNOTE: The GPU has >20%% more load than the CPU. This imbalance causes\n"
849 " performance loss, consider using a shorter cut-off and a finer PME grid.");
855 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
856 (cyc_sum[ewcDOMDEC] > tot*0.1 ||
857 cyc_sum[ewcNS] > tot*0.1))
859 /* Only the sim master calls this function, so always print to stderr */
860 if (wc->wcc[ewcDOMDEC].n == 0)
862 md_print_warn(NULL, fplog,
863 "NOTE: %d %% of the run time was spent in pair search,\n"
864 " you might want to increase nstlist (this has no effect on accuracy)\n",
865 (int)(100*cyc_sum[ewcNS]/tot+0.5));
869 md_print_warn(NULL, fplog,
870 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
871 " %d %% of the run time was spent in pair search,\n"
872 " you might want to increase nstlist (this has no effect on accuracy)\n",
873 (int)(100*cyc_sum[ewcDOMDEC]/tot+0.5),
874 (int)(100*cyc_sum[ewcNS]/tot+0.5));
878 if (cyc_sum[ewcMoveE] > tot*0.05)
880 /* Only the sim master calls this function, so always print to stderr */
881 md_print_warn(NULL, fplog,
882 "NOTE: %d %% of the run time was spent communicating energies,\n"
883 " you might want to use the -gcom option of mdrun\n",
884 (int)(100*cyc_sum[ewcMoveE]/tot+0.5));
888 extern gmx_int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
895 return wc->reset_counters;
898 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_int64_t reset_counters)
905 wc->reset_counters = reset_counters;
908 #ifdef GMX_CYCLE_SUBCOUNTERS
910 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
914 wc->wcsc[ewcs].start = gmx_cycles_read();
918 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
922 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
929 void wallcycle_sub_start(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
932 void wallcycle_sub_stop(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
936 #endif /* GMX_CYCLE_SUBCOUNTERS */