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 "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", "COM pull force",
105 "Write traj.", "Update", "Constraints", "Comm. energies",
106 "Enforced rotation", "Add rot. forces", "Coordinate swapping", "IMD", "Test"
109 static const char *wcsn[ewcsNR] =
111 "DD redist.", "DD NS grid + sort", "DD setup comm.",
112 "DD make top.", "DD make constr.", "DD top. other",
113 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
114 "Bonded F", "Nonbonded F", "Ewald F correction",
115 "NB X buffer ops.", "NB F buffer ops."
118 gmx_bool wallcycle_have_counter(void)
120 return gmx_cycles_have_counter();
123 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec gmx_unused *cr,
124 int nthreads_pp, int nthreads_pme)
129 if (!wallcycle_have_counter())
136 wc->wc_barrier = FALSE;
140 wc->reset_counters = resetstep;
141 wc->nthreads_pp = nthreads_pp;
142 wc->nthreads_pme = nthreads_pme;
143 wc->cycles_sum = NULL;
146 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != NULL)
150 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
152 wc->wc_barrier = TRUE;
153 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
157 snew(wc->wcc, ewcNR);
158 if (getenv("GMX_CYCLE_ALL") != NULL)
162 fprintf(fplog, "\nWill time all the code during the run\n\n");
164 snew(wc->wcc_all, ewcNR*ewcNR);
167 #ifdef GMX_CYCLE_SUBCOUNTERS
168 snew(wc->wcsc, ewcsNR);
178 void wallcycle_destroy(gmx_wallcycle_t wc)
189 if (wc->wcc_all != NULL)
193 #ifdef GMX_CYCLE_SUBCOUNTERS
194 if (wc->wcsc != NULL)
202 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
205 wc->cycle_prev = cycle;
208 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
210 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
211 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
216 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
218 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
220 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
222 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
225 wc->counterlist[wc->count_depth] = ewc;
229 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
233 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
235 if (wc->count_depth < 0)
237 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
239 if (wc->counterlist[wc->count_depth] != ewc)
241 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
242 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
247 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
259 MPI_Barrier(wc->mpi_comm_mygroup);
264 debug_start_check(wc, ewc);
267 cycle = gmx_cycles_read();
268 wc->wcc[ewc].start = cycle;
269 if (wc->wcc_all != NULL)
274 wallcycle_all_start(wc, ewc, cycle);
276 else if (wc->wc_depth == 3)
278 wallcycle_all_stop(wc, ewc, cycle);
283 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
290 wallcycle_start(wc, ewc);
294 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
296 gmx_cycles_t cycle, last;
306 MPI_Barrier(wc->mpi_comm_mygroup);
311 debug_stop_check(wc, ewc);
314 cycle = gmx_cycles_read();
315 last = cycle - wc->wcc[ewc].start;
316 wc->wcc[ewc].c += last;
323 wallcycle_all_stop(wc, ewc, cycle);
325 else if (wc->wc_depth == 2)
327 wallcycle_all_start(wc, ewc, cycle);
334 void wallcycle_reset_all(gmx_wallcycle_t wc)
343 for (i = 0; i < ewcNR; i++)
350 for (i = 0; i < ewcNR*ewcNR; i++)
352 wc->wcc_all[i].n = 0;
353 wc->wcc_all[i].c = 0;
356 #ifdef GMX_CYCLE_SUBCOUNTERS
357 for (i = 0; i < ewcsNR; i++)
365 static gmx_bool is_pme_counter(int ewc)
367 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
370 static gmx_bool is_pme_subcounter(int ewc)
372 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
375 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
378 double cycles[ewcNR+ewcsNR];
379 double cycles_n[ewcNR+ewcsNR], buf[ewcNR+ewcsNR], *cyc_all, *buf_all;
388 snew(wc->cycles_sum, ewcNR+ewcsNR);
392 for (i = 0; i < ewcNR; i++)
394 if (is_pme_counter(i) || (i == ewcRUN && cr->duty == DUTY_PME))
396 wcc[i].c *= wc->nthreads_pme;
400 for (j = 0; j < ewcNR; j++)
402 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pme;
408 wcc[i].c *= wc->nthreads_pp;
412 for (j = 0; j < ewcNR; j++)
414 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pp;
420 if (wcc[ewcDDCOMMLOAD].n > 0)
422 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMLOAD].c;
424 if (wcc[ewcDDCOMMBOUND].n > 0)
426 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMBOUND].c;
428 if (wcc[ewcPME_FFTCOMM].n > 0)
430 wcc[ewcPME_FFT].c -= wcc[ewcPME_FFTCOMM].c;
433 if (cr->npmenodes == 0)
435 /* All nodes do PME (or no PME at all) */
436 if (wcc[ewcPMEMESH].n > 0)
438 wcc[ewcFORCE].c -= wcc[ewcPMEMESH].c;
443 /* The are PME-only nodes */
444 if (wcc[ewcPMEMESH].n > 0)
446 /* This must be a PME only node, calculate the Wait + Comm. time */
447 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
451 /* Store the cycles in a double buffer for summing */
452 for (i = 0; i < ewcNR; i++)
454 cycles_n[i] = (double)wcc[i].n;
455 cycles[i] = (double)wcc[i].c;
458 #ifdef GMX_CYCLE_SUBCOUNTERS
459 for (i = 0; i < ewcsNR; i++)
461 wc->wcsc[i].c *= wc->nthreads_pp;
462 cycles_n[ewcNR+i] = (double)wc->wcsc[i].n;
463 cycles[ewcNR+i] = (double)wc->wcsc[i].c;
471 MPI_Allreduce(cycles_n, buf, nsum, MPI_DOUBLE, MPI_MAX,
473 for (i = 0; i < ewcNR; i++)
475 wcc[i].n = (int)(buf[i] + 0.5);
477 #ifdef GMX_CYCLE_SUBCOUNTERS
478 for (i = 0; i < ewcsNR; i++)
480 wc->wcsc[i].n = (int)(buf[ewcNR+i] + 0.5);
484 MPI_Allreduce(cycles, wc->cycles_sum, nsum, MPI_DOUBLE, MPI_SUM,
487 if (wc->wcc_all != NULL)
489 snew(cyc_all, ewcNR*ewcNR);
490 snew(buf_all, ewcNR*ewcNR);
491 for (i = 0; i < ewcNR*ewcNR; i++)
493 cyc_all[i] = wc->wcc_all[i].c;
495 MPI_Allreduce(cyc_all, buf_all, ewcNR*ewcNR, MPI_DOUBLE, MPI_SUM,
497 for (i = 0; i < ewcNR*ewcNR; i++)
499 wc->wcc_all[i].c = buf_all[i];
508 for (i = 0; i < nsum; i++)
510 wc->cycles_sum[i] = cycles[i];
515 static void print_cycles(FILE *fplog, double c2t, const char *name,
516 int nnodes, int nthreads,
517 int ncalls, double c_sum, double tot)
520 char nthreads_str[6];
528 snprintf(ncalls_str, sizeof(ncalls_str), "%10d", ncalls);
531 snprintf(nnodes_str, sizeof(nnodes_str), "N/A");
535 snprintf(nnodes_str, sizeof(nnodes_str), "%4d", nnodes);
539 snprintf(nthreads_str, sizeof(nthreads_str), "N/A");
543 snprintf(nthreads_str, sizeof(nthreads_str), "%4d", nthreads);
552 /* Convert the cycle count to wallclock time for this task */
555 fprintf(fplog, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
556 name, nnodes_str, nthreads_str, ncalls_str, wallt,
557 c_sum*1e-9, 100*c_sum/tot);
561 static void print_gputimes(FILE *fplog, const char *name,
562 int n, double t, double tot_t)
569 snprintf(num, sizeof(num), "%10d", n);
570 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
575 sprintf(avg_perf, " ");
579 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
580 name, num, t/1000, avg_perf, 100 * t/tot_t);
584 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
585 name, "", t/1000, avg_perf, 100.0);
589 static void print_header(FILE *fplog, int nrank_pp, int nth_pp, int nrank_pme, int nth_pme)
591 int nrank_tot = nrank_pp + nrank_pme;
594 fprintf(fplog, "On %d MPI rank%s", nrank_tot, nrank_tot == 1 ? "" : "s");
597 fprintf(fplog, ", each using %d OpenMP threads", nth_pp);
599 /* Don't report doing PP+PME, because we can't tell here if
600 * this is RF, etc. */
604 fprintf(fplog, "On %d MPI rank%s doing PP", nrank_pp, nrank_pp == 1 ? "" : "s");
607 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pp > 1 ? " each" : "", nth_pp);
609 fprintf(fplog, ", and\non %d MPI rank%s doing PME", nrank_pme, nrank_pme == 1 ? "" : "s");
612 fprintf(fplog, ",%s using %d OpenMP threads", nrank_pme > 1 ? " each" : "", nth_pme);
616 fprintf(fplog, "\n\n");
617 fprintf(fplog, " Computing: Num Num Call Wall time Giga-Cycles\n");
618 fprintf(fplog, " Ranks Threads Count (s) total sum %%\n");
621 void wallcycle_print(FILE *fplog, int nnodes, int npme, double realtime,
622 gmx_wallcycle_t wc, wallclock_gpu_t *gpu_t)
625 double tot, tot_for_pp, tot_for_rest, tot_gpu, tot_cpu_overlap, gpu_cpu_ratio, tot_k;
626 double c2t, c2t_pp, c2t_pme;
627 int i, j, npp, nth_pp, nth_pme, nth_tot;
629 const char *hline = "-----------------------------------------------------------------------------";
636 nth_pp = wc->nthreads_pp;
637 nth_pme = wc->nthreads_pme;
639 cyc_sum = wc->cycles_sum;
642 nth_tot = npp*nth_pp + npme*nth_pme;
644 /* When using PME-only nodes, the next line is valid for both
645 PP-only and PME-only nodes because they started ewcRUN at the
647 tot = cyc_sum[ewcRUN];
650 /* Conversion factor from cycles to seconds */
654 c2t_pp = c2t * nth_tot / (double) (npp*nth_pp);
655 c2t_pme = c2t * nth_tot / (double) (npme*nth_pme);
664 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");
666 print_header(fplog, npp, nth_pp, npme, nth_pme);
668 fprintf(fplog, "%s\n", hline);
669 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
671 if (is_pme_subcounter(i))
673 /* Do not count these at all */
675 else if (npme > 0 && is_pme_counter(i))
677 /* Print timing information for PME-only nodes, but add an
678 * asterisk so the reader of the table can know that the
679 * walltimes are not meant to add up. The asterisk still
680 * fits in the required maximum of 19 characters. */
682 snprintf(buffer, STRLEN, "%s *", wcn[i]);
683 print_cycles(fplog, c2t_pme, buffer,
685 wc->wcc[i].n, cyc_sum[i], tot);
689 /* Print timing information when it is for a PP or PP+PME
691 print_cycles(fplog, c2t_pp, wcn[i],
693 wc->wcc[i].n, cyc_sum[i], tot);
694 tot_for_pp += cyc_sum[i];
697 if (wc->wcc_all != NULL)
699 for (i = 0; i < ewcNR; i++)
701 for (j = 0; j < ewcNR; j++)
703 snprintf(buf, 20, "%-9.9s %-9.9s", wcn[i], wcn[j]);
704 print_cycles(fplog, c2t_pp, buf,
706 wc->wcc_all[i*ewcNR+j].n,
707 wc->wcc_all[i*ewcNR+j].c,
712 tot_for_rest = tot * (npp * nth_pp) / (double) nth_tot;
713 print_cycles(fplog, c2t_pp, "Rest",
715 -1, tot_for_rest - tot_for_pp, tot);
716 fprintf(fplog, "%s\n", hline);
717 print_cycles(fplog, c2t, "Total",
720 fprintf(fplog, "%s\n", hline);
725 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
726 " twice the total reported, but the cycle count total and %% are correct.\n"
730 if (wc->wcc[ewcPMEMESH].n > 0)
732 fprintf(fplog, " Breakdown of PME mesh computation\n");
733 fprintf(fplog, "%s\n", hline);
734 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
736 if (is_pme_subcounter(i))
738 print_cycles(fplog, npme > 0 ? c2t_pme : c2t_pp, wcn[i],
739 npme > 0 ? npme : npp, nth_pme,
740 wc->wcc[i].n, cyc_sum[i], tot);
743 fprintf(fplog, "%s\n", hline);
746 #ifdef GMX_CYCLE_SUBCOUNTERS
747 fprintf(fplog, " Breakdown of PP computation\n");
748 fprintf(fplog, "%s\n", hline);
749 for (i = 0; i < ewcsNR; i++)
751 print_cycles(fplog, c2t_pp, wcsn[i],
753 wc->wcsc[i].n, cyc_sum[ewcNR+i], tot);
755 fprintf(fplog, "%s\n", hline);
758 /* print GPU timing summary */
761 const char *k_log_str[2][2] = {
762 {"Nonbonded F kernel", "Nonbonded F+ene k."},
763 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
766 tot_gpu = gpu_t->pl_h2d_t + gpu_t->nb_h2d_t + gpu_t->nb_d2h_t;
768 /* add up the kernel timings */
770 for (i = 0; i < 2; i++)
772 for (j = 0; j < 2; j++)
774 tot_k += gpu_t->ktime[i][j].t;
779 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
780 if (wc->wcc[ewcPMEMESH].n > 0)
782 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
784 tot_cpu_overlap *= realtime*1000/tot; /* convert s to ms */
786 fprintf(fplog, "\n GPU timings\n%s\n", hline);
787 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
788 fprintf(fplog, "%s\n", hline);
789 print_gputimes(fplog, "Pair list H2D",
790 gpu_t->pl_h2d_c, gpu_t->pl_h2d_t, tot_gpu);
791 print_gputimes(fplog, "X / q H2D",
792 gpu_t->nb_c, gpu_t->nb_h2d_t, tot_gpu);
794 for (i = 0; i < 2; i++)
796 for (j = 0; j < 2; j++)
798 if (gpu_t->ktime[i][j].c)
800 print_gputimes(fplog, k_log_str[i][j],
801 gpu_t->ktime[i][j].c, gpu_t->ktime[i][j].t, tot_gpu);
806 print_gputimes(fplog, "F D2H", gpu_t->nb_c, gpu_t->nb_d2h_t, tot_gpu);
807 fprintf(fplog, "%s\n", hline);
808 print_gputimes(fplog, "Total ", gpu_t->nb_c, tot_gpu, tot_gpu);
809 fprintf(fplog, "%s\n", hline);
811 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
812 fprintf(fplog, "\nForce evaluation time GPU/CPU: %.3f ms/%.3f ms = %.3f\n",
813 tot_gpu/gpu_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
816 /* only print notes related to CPU-GPU load balance with PME */
817 if (wc->wcc[ewcPMEMESH].n > 0)
819 fprintf(fplog, "For optimal performance this ratio should be close to 1!\n");
821 /* print note if the imbalance is high with PME case in which
822 * CPU-GPU load balancing is possible */
823 if (gpu_cpu_ratio < 0.75 || gpu_cpu_ratio > 1.2)
825 /* Only the sim master calls this function, so always print to stderr */
826 if (gpu_cpu_ratio < 0.75)
830 /* The user could have used -notunepme,
831 * but we currently can't check that here.
833 md_print_warn(NULL, fplog,
834 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
835 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
836 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
840 /* We should not end up here, unless the box is
841 * too small for increasing the cut-off for PME tuning.
843 md_print_warn(NULL, fplog,
844 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
845 " performance loss.");
848 if (gpu_cpu_ratio > 1.2)
850 md_print_warn(NULL, fplog,
851 "\nNOTE: The GPU has >20%% more load than the CPU. This imbalance causes\n"
852 " performance loss, consider using a shorter cut-off and a finer PME grid.");
858 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
859 (cyc_sum[ewcDOMDEC] > tot*0.1 ||
860 cyc_sum[ewcNS] > tot*0.1))
862 /* Only the sim master calls this function, so always print to stderr */
863 if (wc->wcc[ewcDOMDEC].n == 0)
865 md_print_warn(NULL, fplog,
866 "NOTE: %d %% of the run time was spent in pair search,\n"
867 " you might want to increase nstlist (this has no effect on accuracy)\n",
868 (int)(100*cyc_sum[ewcNS]/tot+0.5));
872 md_print_warn(NULL, fplog,
873 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
874 " %d %% of the run time was spent in pair search,\n"
875 " you might want to increase nstlist (this has no effect on accuracy)\n",
876 (int)(100*cyc_sum[ewcDOMDEC]/tot+0.5),
877 (int)(100*cyc_sum[ewcNS]/tot+0.5));
881 if (cyc_sum[ewcMoveE] > tot*0.05)
883 /* Only the sim master calls this function, so always print to stderr */
884 md_print_warn(NULL, fplog,
885 "NOTE: %d %% of the run time was spent communicating energies,\n"
886 " you might want to use the -gcom option of mdrun\n",
887 (int)(100*cyc_sum[ewcMoveE]/tot+0.5));
891 extern gmx_int64_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
898 return wc->reset_counters;
901 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_int64_t reset_counters)
908 wc->reset_counters = reset_counters;
911 #ifdef GMX_CYCLE_SUBCOUNTERS
913 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
917 wc->wcsc[ewcs].start = gmx_cycles_read();
921 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
925 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
932 void wallcycle_sub_start(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
935 void wallcycle_sub_stop(gmx_wallcycle_t gmx_unused wc, int gmx_unused ewcs)
939 #endif /* GMX_CYCLE_SUBCOUNTERS */