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"
46 #include "gromacs/legacyheaders/md_logging.h"
47 #include "gromacs/legacyheaders/types/commrec.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/gmxmpi.h"
51 #include "gromacs/utility/smalloc.h"
53 #include "cyclecounter.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 */
69 typedef struct gmx_wallcycle
72 /* variables for testing/debugging */
78 int counterlist[DEPTH_MAX];
82 gmx_cycles_t cycle_prev;
83 gmx_int64_t reset_counters;
85 MPI_Comm mpi_comm_mygroup;
89 #ifdef GMX_CYCLE_SUBCOUNTERS
95 /* Each name should not exceed 19 printing characters
96 (ie. terminating null can be twentieth) */
97 static const char *wcn[ewcNR] =
99 "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
100 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
101 "Comm. coord.", "Born radii", "Force", "Wait + Comm. F", "PME mesh",
102 "PME redist. X/F", "PME spread/gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve LJ", "PME solve Elec",
103 "PME wait for PP", "Wait + Recv. PME F", "Wait GPU nonlocal", "Wait GPU local", "NB X/F buffer ops.",
104 "Vsite spread", "Write traj.", "Update", "Constraints", "Comm. energies",
105 "Enforced rotation", "Add rot. forces", "Coordinate swapping", "IMD", "Test"
108 static const char *wcsn[ewcsNR] =
110 "DD redist.", "DD NS grid + sort", "DD setup comm.",
111 "DD make top.", "DD make constr.", "DD top. other",
112 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
113 "Bonded F", "Nonbonded F", "Ewald F correction",
114 "NB X buffer ops.", "NB F buffer ops."
117 gmx_bool wallcycle_have_counter(void)
119 return gmx_cycles_have_counter();
122 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec gmx_unused *cr,
123 int nthreads_pp, int nthreads_pme)
128 if (!wallcycle_have_counter())
135 wc->wc_barrier = FALSE;
139 wc->reset_counters = resetstep;
140 wc->nthreads_pp = nthreads_pp;
141 wc->nthreads_pme = nthreads_pme;
142 wc->cycles_sum = NULL;
145 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != NULL)
149 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
151 wc->wc_barrier = TRUE;
152 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
156 snew(wc->wcc, ewcNR);
157 if (getenv("GMX_CYCLE_ALL") != NULL)
161 fprintf(fplog, "\nWill time all the code during the run\n\n");
163 snew(wc->wcc_all, ewcNR*ewcNR);
166 #ifdef GMX_CYCLE_SUBCOUNTERS
167 snew(wc->wcsc, ewcsNR);
177 void wallcycle_destroy(gmx_wallcycle_t wc)
188 if (wc->wcc_all != NULL)
192 #ifdef GMX_CYCLE_SUBCOUNTERS
193 if (wc->wcsc != NULL)
201 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
204 wc->cycle_prev = cycle;
207 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
209 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
210 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
215 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
217 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
219 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
221 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
224 wc->counterlist[wc->count_depth] = ewc;
228 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
232 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
234 if (wc->count_depth < 0)
236 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
238 if (wc->counterlist[wc->count_depth] != ewc)
240 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
241 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
246 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
258 MPI_Barrier(wc->mpi_comm_mygroup);
263 debug_start_check(wc, ewc);
266 cycle = gmx_cycles_read();
267 wc->wcc[ewc].start = cycle;
268 if (wc->wcc_all != NULL)
273 wallcycle_all_start(wc, ewc, cycle);
275 else if (wc->wc_depth == 3)
277 wallcycle_all_stop(wc, ewc, cycle);
282 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
289 wallcycle_start(wc, ewc);
293 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
295 gmx_cycles_t cycle, last;
305 MPI_Barrier(wc->mpi_comm_mygroup);
310 debug_stop_check(wc, ewc);
313 cycle = gmx_cycles_read();
314 last = cycle - wc->wcc[ewc].start;
315 wc->wcc[ewc].c += last;
322 wallcycle_all_stop(wc, ewc, cycle);
324 else if (wc->wc_depth == 2)
326 wallcycle_all_start(wc, ewc, cycle);
333 void wallcycle_reset_all(gmx_wallcycle_t wc)
342 for (i = 0; i < ewcNR; i++)
349 for (i = 0; i < ewcNR*ewcNR; i++)
351 wc->wcc_all[i].n = 0;
352 wc->wcc_all[i].c = 0;
355 #ifdef GMX_CYCLE_SUBCOUNTERS
356 for (i = 0; i < ewcsNR; i++)
364 static gmx_bool is_pme_counter(int ewc)
366 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
369 static gmx_bool is_pme_subcounter(int ewc)
371 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
374 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
377 double cycles[ewcNR+ewcsNR];
378 double cycles_n[ewcNR+ewcsNR], buf[ewcNR+ewcsNR], *cyc_all, *buf_all;
387 snew(wc->cycles_sum, ewcNR+ewcsNR);
391 for (i = 0; i < ewcNR; i++)
393 if (is_pme_counter(i) || (i == ewcRUN && cr->duty == DUTY_PME))
395 wcc[i].c *= wc->nthreads_pme;
399 for (j = 0; j < ewcNR; j++)
401 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pme;
407 wcc[i].c *= wc->nthreads_pp;
411 for (j = 0; j < ewcNR; j++)
413 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pp;
419 if (wcc[ewcDDCOMMLOAD].n > 0)
421 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMLOAD].c;
423 if (wcc[ewcDDCOMMBOUND].n > 0)
425 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMBOUND].c;
427 if (wcc[ewcPME_FFTCOMM].n > 0)
429 wcc[ewcPME_FFT].c -= wcc[ewcPME_FFTCOMM].c;
432 if (cr->npmenodes == 0)
434 /* All nodes do PME (or no PME at all) */
435 if (wcc[ewcPMEMESH].n > 0)
437 wcc[ewcFORCE].c -= wcc[ewcPMEMESH].c;
442 /* The are PME-only nodes */
443 if (wcc[ewcPMEMESH].n > 0)
445 /* This must be a PME only node, calculate the Wait + Comm. time */
446 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
450 /* Store the cycles in a double buffer for summing */
451 for (i = 0; i < ewcNR; i++)
453 cycles_n[i] = (double)wcc[i].n;
454 cycles[i] = (double)wcc[i].c;
457 #ifdef GMX_CYCLE_SUBCOUNTERS
458 for (i = 0; i < ewcsNR; i++)
460 wc->wcsc[i].c *= wc->nthreads_pp;
461 cycles_n[ewcNR+i] = (double)wc->wcsc[i].n;
462 cycles[ewcNR+i] = (double)wc->wcsc[i].c;
470 MPI_Allreduce(cycles_n, buf, nsum, MPI_DOUBLE, MPI_MAX,
472 for (i = 0; i < ewcNR; i++)
474 wcc[i].n = (int)(buf[i] + 0.5);
476 #ifdef GMX_CYCLE_SUBCOUNTERS
477 for (i = 0; i < ewcsNR; i++)
479 wc->wcsc[i].n = (int)(buf[ewcNR+i] + 0.5);
483 MPI_Allreduce(cycles, wc->cycles_sum, nsum, MPI_DOUBLE, MPI_SUM,
486 if (wc->wcc_all != NULL)
488 snew(cyc_all, ewcNR*ewcNR);
489 snew(buf_all, ewcNR*ewcNR);
490 for (i = 0; i < ewcNR*ewcNR; i++)
492 cyc_all[i] = wc->wcc_all[i].c;
494 MPI_Allreduce(cyc_all, buf_all, ewcNR*ewcNR, MPI_DOUBLE, MPI_SUM,
496 for (i = 0; i < ewcNR*ewcNR; i++)
498 wc->wcc_all[i].c = buf_all[i];
507 for (i = 0; i < nsum; i++)
509 wc->cycles_sum[i] = cycles[i];
514 static void print_cycles(FILE *fplog, double c2t, const char *name,
515 int nnodes, int nthreads,
516 int ncalls, double c_sum, double tot)
519 char nthreads_str[6];
527 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
530 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
534 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
538 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
542 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
551 /* Convert the cycle count to wallclock time for this task */
554 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
555 name, nnodes_str, nthreads_str, ncalls_str, wallt,
556 c_sum*1e-9, 100*c_sum/tot);
560 static void print_gputimes(FILE *fplog, const char *name,
561 int n, double t, double tot_t)
568 snprintf(num, sizeof(num), "%10d", n);
569 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
574 sprintf(avg_perf, " ");
578 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
579 name, num, t/1000, avg_perf, 100 * t/tot_t);
583 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
584 name, "", t/1000, avg_perf, 100.0);
588 static void print_header(FILE *fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
590 int nrank_tot = nrank_pp + nrank_pme;
593 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
596 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
598 /* Don't report doing PP+PME, because we can't tell here if
599 * this is RF, etc. */
603 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
606 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
608 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
611 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
615 fprintf(fplog, "\n\n");
616 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
617 fprintf(fplog, " Nodes Threads Count (s) total sum %%\n");
620 void wallcycle_print(FILE *fplog, int nnodes, int npme, double realtime,
621 gmx_wallcycle_t wc, wallclock_gpu_t *gpu_t)
624 double tot, tot_for_pp, tot_for_rest, tot_gpu, tot_cpu_overlap, gpu_cpu_ratio, tot_k;
625 double c2t, c2t_pp, c2t_pme;
626 int i, j, npp, nth_pp, nth_pme, nth_tot;
628 const char *hline = "-----------------------------------------------------------------------------";
635 nth_pp = wc->nthreads_pp;
636 nth_pme = wc->nthreads_pme;
638 cyc_sum = wc->cycles_sum;
641 nth_tot = npp*nth_pp + npme*nth_pme;
643 /* When using PME-only nodes, the next line is valid for both
644 PP-only and PME-only nodes because they started ewcRUN at the
646 tot = cyc_sum[ewcRUN];
649 /* Conversion factor from cycles to seconds */
653 c2t_pp = c2t * nth_tot / (double) (npp*nth_pp);
654 c2t_pme = c2t * nth_tot / (double) (npme*nth_pme);
663 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");
665 print_header(fplog, npp, nth_pp, npme, nth_pme);
667 fprintf(fplog, "%s\n", hline);
668 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
670 if (is_pme_subcounter(i))
672 /* Do not count these at all */
674 else if (npme > 0 && is_pme_counter(i))
676 /* Print timing information for PME-only nodes, but add an
677 * asterisk so the reader of the table can know that the
678 * walltimes are not meant to add up. The asterisk still
679 * fits in the required maximum of 19 characters. */
681 snprintf(buffer, STRLEN, "%s *", wcn[i]);
682 print_cycles(fplog, c2t_pme, buffer,
684 wc->wcc[i].n, cyc_sum[i], tot);
688 /* Print timing information when it is for a PP or PP+PME
690 print_cycles(fplog, c2t_pp, wcn[i],
692 wc->wcc[i].n, cyc_sum[i], tot);
693 tot_for_pp += cyc_sum[i];
696 if (wc->wcc_all != NULL)
698 for (i = 0; i < ewcNR; i++)
700 for (j = 0; j < ewcNR; j++)
702 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
703 print_cycles(fplog, c2t_pp, buf,
705 wc->wcc_all[i*ewcNR+j].n,
706 wc->wcc_all[i*ewcNR+j].c,
711 tot_for_rest = tot * (npp * nth_pp) / (double) nth_tot;
712 print_cycles(fplog, c2t_pp, "Rest",
714 -1, tot_for_rest - tot_for_pp, tot);
715 fprintf(fplog, "%s\n", hline);
716 print_cycles(fplog, c2t, "Total",
719 fprintf(fplog, "%s\n", hline);
724 "(*) Note that with separate PME nodes, the walltime column actually sums to\n"
725 " twice the total reported, but the cycle count total and %% are correct.\n"
729 if (wc->wcc[ewcPMEMESH].n > 0)
731 fprintf(fplog, " Breakdown of PME mesh computation\n");
732 fprintf(fplog, "%s\n", hline);
733 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
735 if (is_pme_subcounter(i))
737 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i],
738 npme > 0 ? npme : npp, nth_pme,
739 wc->wcc[i].n, cyc_sum[i], tot);
742 fprintf(fplog, "%s\n", hline);
745 #ifdef GMX_CYCLE_SUBCOUNTERS
746 fprintf(fplog, " Breakdown of PP computation\n");
747 fprintf(fplog, "%s\n", hline);
748 for (i = 0; i < ewcsNR; i++)
750 print_cycles(fplog, c2t_pp, wcsn[i],
752 wc->wcsc[i].n, cyc_sum[ewcNR+i], tot);
754 fprintf(fplog, "%s\n", hline);
757 /* print GPU timing summary */
760 const char *k_log_str[2][2] = {
761 {"Nonbonded F kernel", "Nonbonded F+ene k."},
762 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
765 tot_gpu = gpu_t->pl_h2d_t + gpu_t->nb_h2d_t + gpu_t->nb_d2h_t;
767 /* add up the kernel timings */
769 for (i = 0; i < 2; i++)
771 for (j = 0; j < 2; j++)
773 tot_k += gpu_t->ktime[i][j].t;
778 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
779 if (wc->wcc[ewcPMEMESH].n > 0)
781 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
783 tot_cpu_overlap *= realtime*1000/tot; /* convert s to ms */
785 fprintf(fplog, "\n GPU timings\n%s\n", hline);
786 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
787 fprintf(fplog, "%s\n", hline);
788 print_gputimes(fplog, "Pair list H2D",
789 gpu_t->pl_h2d_c, gpu_t->pl_h2d_t, tot_gpu);
790 print_gputimes(fplog, "X / q H2D",
791 gpu_t->nb_c, gpu_t->nb_h2d_t, tot_gpu);
793 for (i = 0; i < 2; i++)
795 for (j = 0; j < 2; j++)
797 if (gpu_t->ktime[i][j].c)
799 print_gputimes(fplog, k_log_str[i][j],
800 gpu_t->ktime[i][j].c, gpu_t->ktime[i][j].t, tot_gpu);
805 print_gputimes(fplog, "F D2H", gpu_t->nb_c, gpu_t->nb_d2h_t, tot_gpu);
806 fprintf(fplog, "%s\n", hline);
807 print_gputimes(fplog, "Total ", gpu_t->nb_c, tot_gpu, tot_gpu);
808 fprintf(fplog, "%s\n", hline);
810 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
811 fprintf(fplog, "\nForce evaluation time GPU/CPU: %.3f ms/%.3f ms = %.3f\n",
812 tot_gpu/gpu_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
815 /* only print notes related to CPU-GPU load balance with PME */
816 if (wc->wcc[ewcPMEMESH].n > 0)
818 fprintf(fplog, "For optimal performance this ratio should be close to 1!\n");
820 /* print note if the imbalance is high with PME case in which
821 * CPU-GPU load balancing is possible */
822 if (gpu_cpu_ratio < 0.75 || gpu_cpu_ratio > 1.2)
824 /* Only the sim master calls this function, so always print to stderr */
825 if (gpu_cpu_ratio < 0.75)
829 /* The user could have used -notunepme,
830 * but we currently can't check that here.
832 md_print_warn(NULL, fplog,
833 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
834 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
835 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
839 /* We should not end up here, unless the box is
840 * too small for increasing the cut-off for PME tuning.
842 md_print_warn(NULL, fplog,
843 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
844 " performance loss.");
847 if (gpu_cpu_ratio > 1.2)
849 md_print_warn(NULL, fplog,
850 "\nNOTE: The GPU has >20%% more load than the CPU. This imbalance causes\n"
851 " performance loss, consider using a shorter cut-off and a finer PME grid.");
857 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
858 (cyc_sum[ewcDOMDEC] > tot*0.1 ||
859 cyc_sum[ewcNS] > tot*0.1))
861 /* Only the sim master calls this function, so always print to stderr */
862 if (wc->wcc[ewcDOMDEC].n == 0)
864 md_print_warn(NULL, fplog,
865 "NOTE: %d %% of the run time was spent in pair search,\n"
866 " you might want to increase nstlist (this has no effect on accuracy)\n",
867 (int)(100*cyc_sum[ewcNS]/tot+0.5));
871 md_print_warn(NULL, fplog,
872 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
873 " %d %% of the run time was spent in pair search,\n"
874 " you might want to increase nstlist (this has no effect on accuracy)\n",
875 (int)(100*cyc_sum[ewcDOMDEC]/tot+0.5),
876 (int)(100*cyc_sum[ewcNS]/tot+0.5));
880 if (cyc_sum[ewcMoveE] > tot*0.05)
882 /* Only the sim master calls this function, so always print to stderr */
883 md_print_warn(NULL, fplog,
884 "NOTE: %d %% of the run time was spent communicating energies,\n"
885 " you might want to use the -gcom option of mdrun\n",
886 (int)(100*cyc_sum[ewcMoveE]/tot+0.5));
890 extern gmx_int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
897 return wc->reset_counters;
900 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_int64_t reset_counters)
907 wc->reset_counters = reset_counters;
910 #ifdef GMX_CYCLE_SUBCOUNTERS
912 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
916 wc->wcsc[ewcs].start = gmx_cycles_read();
920 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
924 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
931 void wallcycle_sub_start(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
934 void wallcycle_sub_stop(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
938 #endif /* GMX_CYCLE_SUBCOUNTERS */