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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gmx_wallcycle.h"
46 #include "gmx_cyclecounter.h"
48 #include "gmx_fatal.h"
49 #include "md_logging.h"
59 /* DEBUG_WCYCLE adds consistency checking for the counters.
60 * It checks if you stop a counter different from the last
61 * one that was opened and if you do nest too deep.
63 /* #define DEBUG_WCYCLE */
73 typedef struct gmx_wallcycle
76 /* variables for testing/debugging */
82 int counterlist[DEPTH_MAX];
86 gmx_cycles_t cycle_prev;
87 gmx_large_int_t reset_counters;
89 MPI_Comm mpi_comm_mygroup;
93 #ifdef GMX_CYCLE_SUBCOUNTERS
99 /* Each name should not exceed 19 characters */
100 static const char *wcn[ewcNR] =
102 "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
103 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
104 "Comm. coord.", "Born radii", "Force", "Wait + Comm. F", "PME mesh",
105 "PME redist. X/F", "PME spread/gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve",
106 "PME wait for PP", "Wait + Recv. PME F", "Wait GPU nonlocal", "Wait GPU local", "NB X/F buffer ops.",
107 "Vsite spread", "Write traj.", "Update", "Constraints", "Comm. energies",
108 "Enforced rotation", "Add rot. forces", "Test"
111 static const char *wcsn[ewcsNR] =
113 "DD redist.", "DD NS grid + sort", "DD setup comm.",
114 "DD make top.", "DD make constr.", "DD top. other",
115 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
116 "Bonded F", "Nonbonded F", "Ewald F correction",
117 "NB X buffer ops.", "NB F buffer ops."
120 gmx_bool wallcycle_have_counter(void)
122 return gmx_cycles_have_counter();
125 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec *cr,
126 int nthreads_pp, int nthreads_pme)
131 if (!wallcycle_have_counter())
138 wc->wc_barrier = FALSE;
142 wc->reset_counters = resetstep;
143 wc->nthreads_pp = nthreads_pp;
144 wc->nthreads_pme = nthreads_pme;
145 wc->cycles_sum = NULL;
148 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != NULL)
152 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
154 wc->wc_barrier = TRUE;
155 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
159 snew(wc->wcc, ewcNR);
160 if (getenv("GMX_CYCLE_ALL") != NULL)
164 fprintf(fplog, "\nWill time all the code during the run\n\n");
166 snew(wc->wcc_all, ewcNR*ewcNR);
169 #ifdef GMX_CYCLE_SUBCOUNTERS
170 snew(wc->wcsc, ewcsNR);
180 void wallcycle_destroy(gmx_wallcycle_t wc)
191 if (wc->wcc_all != NULL)
195 #ifdef GMX_CYCLE_SUBCOUNTERS
196 if (wc->wcsc != NULL)
204 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
207 wc->cycle_prev = cycle;
210 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
212 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
213 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
218 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
220 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
222 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
224 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
227 wc->counterlist[wc->count_depth] = ewc;
231 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
235 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
237 if (wc->count_depth < 0)
239 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
241 if (wc->counterlist[wc->count_depth] != ewc)
243 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
244 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
249 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
261 MPI_Barrier(wc->mpi_comm_mygroup);
266 debug_start_check(wc, ewc);
269 cycle = gmx_cycles_read();
270 wc->wcc[ewc].start = cycle;
271 if (wc->wcc_all != NULL)
276 wallcycle_all_start(wc, ewc, cycle);
278 else if (wc->wc_depth == 3)
280 wallcycle_all_stop(wc, ewc, cycle);
285 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
292 wallcycle_start(wc, ewc);
296 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
298 gmx_cycles_t cycle, last;
308 MPI_Barrier(wc->mpi_comm_mygroup);
313 debug_stop_check(wc, ewc);
316 cycle = gmx_cycles_read();
317 last = cycle - wc->wcc[ewc].start;
318 wc->wcc[ewc].c += last;
325 wallcycle_all_stop(wc, ewc, cycle);
327 else if (wc->wc_depth == 2)
329 wallcycle_all_start(wc, ewc, cycle);
336 void wallcycle_reset_all(gmx_wallcycle_t wc)
345 for (i = 0; i < ewcNR; i++)
352 for (i = 0; i < ewcNR*ewcNR; i++)
354 wc->wcc_all[i].n = 0;
355 wc->wcc_all[i].c = 0;
358 #ifdef GMX_CYCLE_SUBCOUNTERS
359 for (i = 0; i < ewcsNR; i++)
367 static gmx_bool is_pme_counter(int ewc)
369 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
372 static gmx_bool is_pme_subcounter(int ewc)
374 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
377 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
381 double cycles_n[ewcNR+ewcsNR], buf[ewcNR+ewcsNR], *cyc_all, *buf_all;
390 snew(wc->cycles_sum, ewcNR+ewcsNR);
391 cycles = wc->cycles_sum;
395 for (i = 0; i < ewcNR; i++)
397 if (is_pme_counter(i) || (i == ewcRUN && cr->duty == DUTY_PME))
399 wcc[i].c *= wc->nthreads_pme;
403 for (j = 0; j < ewcNR; j++)
405 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pme;
411 wcc[i].c *= wc->nthreads_pp;
415 for (j = 0; j < ewcNR; j++)
417 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pp;
423 if (wcc[ewcDDCOMMLOAD].n > 0)
425 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMLOAD].c;
427 if (wcc[ewcDDCOMMBOUND].n > 0)
429 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMBOUND].c;
431 if (wcc[ewcPME_FFTCOMM].n > 0)
433 wcc[ewcPME_FFT].c -= wcc[ewcPME_FFTCOMM].c;
436 if (cr->npmenodes == 0)
438 /* All nodes do PME (or no PME at all) */
439 if (wcc[ewcPMEMESH].n > 0)
441 wcc[ewcFORCE].c -= wcc[ewcPMEMESH].c;
446 /* The are PME-only nodes */
447 if (wcc[ewcPMEMESH].n > 0)
449 /* This must be a PME only node, calculate the Wait + Comm. time */
450 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
454 /* Store the cycles in a double buffer for summing */
455 for (i = 0; i < ewcNR; i++)
457 cycles_n[i] = (double)wcc[i].n;
458 cycles[i] = (double)wcc[i].c;
461 #ifdef GMX_CYCLE_SUBCOUNTERS
462 for (i = 0; i < ewcsNR; i++)
464 wc->wcsc[i].c *= wc->nthreads_pp;
465 cycles_n[ewcNR+i] = (double)wc->wcsc[i].n;
466 cycles[ewcNR+i] = (double)wc->wcsc[i].c;
474 MPI_Allreduce(cycles_n, buf, nsum, MPI_DOUBLE, MPI_MAX,
476 for (i = 0; i < ewcNR; i++)
478 wcc[i].n = (int)(buf[i] + 0.5);
480 #ifdef GMX_CYCLE_SUBCOUNTERS
481 for (i = 0; i < ewcsNR; i++)
483 wc->wcsc[i].n = (int)(buf[ewcNR+i] + 0.5);
487 MPI_Allreduce(cycles, buf, nsum, MPI_DOUBLE, MPI_SUM,
489 for (i = 0; i < nsum; i++)
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 static void print_cycles(FILE *fplog, double c2t, const char *name,
517 int nnodes, int nthreads,
518 int n, double c, double tot)
528 snprintf(num, sizeof(num), "%10d", n);
531 snprintf(thstr, sizeof(thstr), "N/A");
535 snprintf(thstr, sizeof(thstr), "%4d", nthreads);
543 /* Convert the cycle count to wallclock time for this task */
546 /* Cycle count has been multiplied by the thread count,
547 * correct for the number of threads used.
549 wallt = c*c2t*nthreads_tot/(double)(nnodes*nthreads);
553 /* nthreads=-1 signals total run time, no correction required */
556 fprintf(fplog, " %-19s %4d %4s %10s %10.3f %12.3f %5.1f\n",
557 name, nnodes, thstr, num, wallt, c*1e-9, 100*c/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 void wallcycle_print(FILE *fplog, int nnodes, int npme, double realtime,
590 gmx_wallcycle_t wc, wallclock_gpu_t *gpu_t)
593 double c2t, tot, tot_gpu, tot_cpu_overlap, gpu_cpu_ratio, sum, tot_k;
594 int i, j, npp, nth_pp, nth_pme, nth_tot;
596 const char *hline = "-----------------------------------------------------------------------------";
603 nth_pp = wc->nthreads_pp;
604 nth_pme = wc->nthreads_pme;
606 cycles = wc->cycles_sum;
612 nth_tot = npp*nth_pp + npme*nth_pme;
619 nth_tot = npp*nth_pp;
622 tot = cycles[ewcRUN];
624 /* Conversion factor from cycles to seconds */
634 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");
636 fprintf(fplog, " Computing: Nodes Th. Count Wall t (s) G-Cycles %c\n", '%');
637 fprintf(fplog, "%s\n", hline);
639 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
641 if (!is_pme_subcounter(i))
643 print_cycles(fplog, c2t, wcn[i], nth_tot,
644 is_pme_counter(i) ? npme : npp,
645 is_pme_counter(i) ? nth_pme : nth_pp,
646 wc->wcc[i].n, cycles[i], tot);
650 if (wc->wcc_all != NULL)
652 for (i = 0; i < ewcNR; i++)
654 for (j = 0; j < ewcNR; j++)
656 snprintf(buf, 9, "%-9s", wcn[i]);
658 snprintf(buf+10, 9, "%-9s", wcn[j]);
660 print_cycles(fplog, c2t, buf, nth_tot,
661 is_pme_counter(i) ? npme : npp,
662 is_pme_counter(i) ? nth_pme : nth_pp,
663 wc->wcc_all[i*ewcNR+j].n,
664 wc->wcc_all[i*ewcNR+j].c,
669 print_cycles(fplog, c2t, "Rest", nth_tot, npp, -1, 0, tot-sum, tot);
670 fprintf(fplog, "%s\n", hline);
671 print_cycles(fplog, c2t, "Total", nth_tot, nnodes, -1, 0, tot, tot);
672 fprintf(fplog, "%s\n", hline);
674 if (wc->wcc[ewcPMEMESH].n > 0)
676 fprintf(fplog, "%s\n", hline);
677 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
679 if (is_pme_subcounter(i))
681 print_cycles(fplog, c2t, wcn[i], nth_tot,
682 is_pme_counter(i) ? npme : npp,
683 is_pme_counter(i) ? nth_pme : nth_pp,
684 wc->wcc[i].n, cycles[i], tot);
687 fprintf(fplog, "%s\n", hline);
690 #ifdef GMX_CYCLE_SUBCOUNTERS
691 fprintf(fplog, "%s\n", hline);
692 for (i = 0; i < ewcsNR; i++)
694 print_cycles(fplog, c2t, wcsn[i], nth_tot, npp, nth_pp,
695 wc->wcsc[i].n, cycles[ewcNR+i], tot);
697 fprintf(fplog, "%s\n", hline);
700 /* print GPU timing summary */
703 const char *k_log_str[2][2] = {
704 {"Nonbonded F kernel", "Nonbonded F+ene k."},
705 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
708 tot_gpu = gpu_t->pl_h2d_t + gpu_t->nb_h2d_t + gpu_t->nb_d2h_t;
710 /* add up the kernel timings */
712 for (i = 0; i < 2; i++)
714 for (j = 0; j < 2; j++)
716 tot_k += gpu_t->ktime[i][j].t;
721 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
722 if (wc->wcc[ewcPMEMESH].n > 0)
724 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
726 tot_cpu_overlap *= c2t * 1000; /* convert s to ms */
728 fprintf(fplog, "\n GPU timings\n%s\n", hline);
729 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
730 fprintf(fplog, "%s\n", hline);
731 print_gputimes(fplog, "Pair list H2D",
732 gpu_t->pl_h2d_c, gpu_t->pl_h2d_t, tot_gpu);
733 print_gputimes(fplog, "X / q H2D",
734 gpu_t->nb_c, gpu_t->nb_h2d_t, tot_gpu);
736 for (i = 0; i < 2; i++)
738 for (j = 0; j < 2; j++)
740 if (gpu_t->ktime[i][j].c)
742 print_gputimes(fplog, k_log_str[i][j],
743 gpu_t->ktime[i][j].c, gpu_t->ktime[i][j].t, tot_gpu);
748 print_gputimes(fplog, "F D2H", gpu_t->nb_c, gpu_t->nb_d2h_t, tot_gpu);
749 fprintf(fplog, "%s\n", hline);
750 print_gputimes(fplog, "Total ", gpu_t->nb_c, tot_gpu, tot_gpu);
751 fprintf(fplog, "%s\n", hline);
753 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
754 fprintf(fplog, "\nForce evaluation time GPU/CPU: %.3f ms/%.3f ms = %.3f\n",
755 tot_gpu/gpu_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
758 /* only print notes related to CPU-GPU load balance with PME */
759 if (wc->wcc[ewcPMEMESH].n > 0)
761 fprintf(fplog, "For optimal performance this ratio should be close to 1!\n");
763 /* print note if the imbalance is high with PME case in which
764 * CPU-GPU load balancing is possible */
765 if (gpu_cpu_ratio < 0.75 || gpu_cpu_ratio > 1.2)
767 /* Only the sim master calls this function, so always print to stderr */
768 if (gpu_cpu_ratio < 0.75)
772 /* The user could have used -notunepme,
773 * but we currently can't check that here.
775 md_print_warn(NULL, fplog,
776 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
777 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
778 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
782 /* We should not end up here, unless the box is
783 * too small for increasing the cut-off for PME tuning.
785 md_print_warn(NULL, fplog,
786 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
787 " performance loss.");
790 if (gpu_cpu_ratio > 1.2)
792 md_print_warn(NULL, fplog,
793 "\nNOTE: The GPU has >20%% more load than the CPU. This imbalance causes\n"
794 " performance loss, consider using a shorter cut-off and a finer PME grid.");
800 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
801 (cycles[ewcDOMDEC] > tot*0.1 ||
802 cycles[ewcNS] > tot*0.1))
804 /* Only the sim master calls this function, so always print to stderr */
805 if (wc->wcc[ewcDOMDEC].n == 0)
807 md_print_warn(NULL, fplog,
808 "NOTE: %d %% of the run time was spent in pair search,\n"
809 " you might want to increase nstlist (this has no effect on accuracy)\n",
810 (int)(100*cycles[ewcNS]/tot+0.5));
814 md_print_warn(NULL, fplog,
815 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
816 " %d %% of the run time was spent in pair search,\n"
817 " you might want to increase nstlist (this has no effect on accuracy)\n",
818 (int)(100*cycles[ewcDOMDEC]/tot+0.5),
819 (int)(100*cycles[ewcNS]/tot+0.5));
823 if (cycles[ewcMoveE] > tot*0.05)
825 /* Only the sim master calls this function, so always print to stderr */
826 md_print_warn(NULL, fplog,
827 "NOTE: %d %% of the run time was spent communicating energies,\n"
828 " you might want to use the -gcom option of mdrun\n",
829 (int)(100*cycles[ewcMoveE]/tot+0.5));
833 extern gmx_large_int_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
840 return wc->reset_counters;
843 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_large_int_t reset_counters)
850 wc->reset_counters = reset_counters;
853 #ifdef GMX_CYCLE_SUBCOUNTERS
855 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
859 wc->wcsc[ewcs].start = gmx_cycles_read();
863 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
867 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
872 #endif /* GMX_CYCLE_SUBCOUNTERS */