1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2008, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
42 #include "gmx_wallcycle.h"
43 #include "gmx_cyclecounter.h"
45 #include "gmx_fatal.h"
46 #include "md_logging.h"
49 #include "gromacs/utility/gmxmpi.h"
51 /* DEBUG_WCYCLE adds consistency checking for the counters.
52 * It checks if you stop a counter different from the last
53 * one that was opened and if you do nest too deep.
55 /* #define DEBUG_WCYCLE */
65 typedef struct gmx_wallcycle
68 /* variables for testing/debugging */
74 int counterlist[DEPTH_MAX];
78 gmx_cycles_t cycle_prev;
79 gmx_large_int_t reset_counters;
81 MPI_Comm mpi_comm_mygroup;
85 #ifdef GMX_CYCLE_SUBCOUNTERS
91 /* Each name should not exceed 19 characters */
92 static const char *wcn[ewcNR] =
94 "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
95 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
96 "Comm. coord.", "Born radii", "Force", "Wait + Comm. F", "PME mesh",
97 "PME redist. X/F", "PME spread/gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve",
98 "PME wait for PP", "Wait + Recv. PME F", "Wait GPU nonlocal", "Wait GPU local", "NB X/F buffer ops.",
99 "Vsite spread", "Write traj.", "Update", "Constraints", "Comm. energies",
100 "Enforced rotation", "Add rot. forces", "Test"
103 static const char *wcsn[ewcsNR] =
105 "DD redist.", "DD NS grid + sort", "DD setup comm.",
106 "DD make top.", "DD make constr.", "DD top. other",
107 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
108 "Bonded F", "Nonbonded F", "Ewald F correction",
109 "NB X buffer ops.", "NB F buffer ops."
112 gmx_bool wallcycle_have_counter(void)
114 return gmx_cycles_have_counter();
117 gmx_wallcycle_t wallcycle_init(FILE *fplog, int resetstep, t_commrec *cr,
118 int nthreads_pp, int nthreads_pme)
123 if (!wallcycle_have_counter())
130 wc->wc_barrier = FALSE;
134 wc->reset_counters = resetstep;
135 wc->nthreads_pp = nthreads_pp;
136 wc->nthreads_pme = nthreads_pme;
137 wc->cycles_sum = NULL;
140 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != NULL)
144 fprintf(fplog, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
146 wc->wc_barrier = TRUE;
147 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
151 snew(wc->wcc, ewcNR);
152 if (getenv("GMX_CYCLE_ALL") != NULL)
156 fprintf(fplog, "\nWill time all the code during the run\n\n");
158 snew(wc->wcc_all, ewcNR*ewcNR);
161 #ifdef GMX_CYCLE_SUBCOUNTERS
162 snew(wc->wcsc, ewcsNR);
172 void wallcycle_destroy(gmx_wallcycle_t wc)
183 if (wc->wcc_all != NULL)
187 #ifdef GMX_CYCLE_SUBCOUNTERS
188 if (wc->wcsc != NULL)
196 static void wallcycle_all_start(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
199 wc->cycle_prev = cycle;
202 static void wallcycle_all_stop(gmx_wallcycle_t wc, int ewc, gmx_cycles_t cycle)
204 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
205 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
210 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
212 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
214 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
216 gmx_fatal(FARGS, "wallcycle counter depth out of range: %d",
219 wc->counterlist[wc->count_depth] = ewc;
223 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
227 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
229 if (wc->count_depth < 0)
231 gmx_fatal(FARGS, "wallcycle counter depth out of range when stopping %s: %d", wcn[ewc], wc->count_depth);
233 if (wc->counterlist[wc->count_depth] != ewc)
235 gmx_fatal(FARGS, "wallcycle mismatch at stop, start %s, stop %s",
236 wcn[wc->counterlist[wc->count_depth]], wcn[ewc]);
241 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
253 MPI_Barrier(wc->mpi_comm_mygroup);
258 debug_start_check(wc, ewc);
261 cycle = gmx_cycles_read();
262 wc->wcc[ewc].start = cycle;
263 if (wc->wcc_all != NULL)
268 wallcycle_all_start(wc, ewc, cycle);
270 else if (wc->wc_depth == 3)
272 wallcycle_all_stop(wc, ewc, cycle);
277 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
284 wallcycle_start(wc, ewc);
288 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
290 gmx_cycles_t cycle, last;
300 MPI_Barrier(wc->mpi_comm_mygroup);
305 debug_stop_check(wc, ewc);
308 cycle = gmx_cycles_read();
309 last = cycle - wc->wcc[ewc].start;
310 wc->wcc[ewc].c += last;
317 wallcycle_all_stop(wc, ewc, cycle);
319 else if (wc->wc_depth == 2)
321 wallcycle_all_start(wc, ewc, cycle);
328 void wallcycle_reset_all(gmx_wallcycle_t wc)
337 for (i = 0; i < ewcNR; i++)
344 for (i = 0; i < ewcNR*ewcNR; i++)
346 wc->wcc_all[i].n = 0;
347 wc->wcc_all[i].c = 0;
350 #ifdef GMX_CYCLE_SUBCOUNTERS
351 for (i = 0; i < ewcsNR; i++)
359 static gmx_bool is_pme_counter(int ewc)
361 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
364 static gmx_bool is_pme_subcounter(int ewc)
366 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
369 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
373 double cycles_n[ewcNR+ewcsNR], buf[ewcNR+ewcsNR], *cyc_all, *buf_all;
382 snew(wc->cycles_sum, ewcNR+ewcsNR);
383 cycles = wc->cycles_sum;
387 for (i = 0; i < ewcNR; i++)
389 if (is_pme_counter(i) || (i == ewcRUN && cr->duty == DUTY_PME))
391 wcc[i].c *= wc->nthreads_pme;
395 for (j = 0; j < ewcNR; j++)
397 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pme;
403 wcc[i].c *= wc->nthreads_pp;
407 for (j = 0; j < ewcNR; j++)
409 wc->wcc_all[i*ewcNR+j].c *= wc->nthreads_pp;
415 if (wcc[ewcDDCOMMLOAD].n > 0)
417 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMLOAD].c;
419 if (wcc[ewcDDCOMMBOUND].n > 0)
421 wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMBOUND].c;
423 if (wcc[ewcPME_FFTCOMM].n > 0)
425 wcc[ewcPME_FFT].c -= wcc[ewcPME_FFTCOMM].c;
428 if (cr->npmenodes == 0)
430 /* All nodes do PME (or no PME at all) */
431 if (wcc[ewcPMEMESH].n > 0)
433 wcc[ewcFORCE].c -= wcc[ewcPMEMESH].c;
438 /* The are PME-only nodes */
439 if (wcc[ewcPMEMESH].n > 0)
441 /* This must be a PME only node, calculate the Wait + Comm. time */
442 wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
446 /* Store the cycles in a double buffer for summing */
447 for (i = 0; i < ewcNR; i++)
449 cycles_n[i] = (double)wcc[i].n;
450 cycles[i] = (double)wcc[i].c;
453 #ifdef GMX_CYCLE_SUBCOUNTERS
454 for (i = 0; i < ewcsNR; i++)
456 wc->wcsc[i].c *= wc->nthreads_pp;
457 cycles_n[ewcNR+i] = (double)wc->wcsc[i].n;
458 cycles[ewcNR+i] = (double)wc->wcsc[i].c;
466 MPI_Allreduce(cycles_n, buf, nsum, MPI_DOUBLE, MPI_MAX,
468 for (i = 0; i < ewcNR; i++)
470 wcc[i].n = (int)(buf[i] + 0.5);
472 #ifdef GMX_CYCLE_SUBCOUNTERS
473 for (i = 0; i < ewcsNR; i++)
475 wc->wcsc[i].n = (int)(buf[ewcNR+i] + 0.5);
479 MPI_Allreduce(cycles, buf, nsum, MPI_DOUBLE, MPI_SUM,
481 for (i = 0; i < nsum; i++)
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 static void print_cycles(FILE *fplog, double c2t, const char *name,
509 int nnodes, int nthreads,
510 int n, double c, double tot)
520 snprintf(num, sizeof(num), "%10d", n);
523 snprintf(thstr, sizeof(thstr), "N/A");
527 snprintf(thstr, sizeof(thstr), "%4d", nthreads);
535 /* Convert the cycle count to wallclock time for this task */
538 /* Cycle count has been multiplied by the thread count,
539 * correct for the number of threads used.
541 wallt = c*c2t*nthreads_tot/(double)(nnodes*nthreads);
545 /* nthreads=-1 signals total run time, no correction required */
548 fprintf(fplog, " %-19s %4d %4s %10s %10.3f %12.3f %5.1f\n",
549 name, nnodes, thstr, num, wallt, c*1e-9, 100*c/tot);
553 static void print_gputimes(FILE *fplog, const char *name,
554 int n, double t, double tot_t)
561 snprintf(num, sizeof(num), "%10d", n);
562 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
567 sprintf(avg_perf, " ");
571 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
572 name, num, t/1000, avg_perf, 100 * t/tot_t);
576 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
577 name, "", t/1000, avg_perf, 100.0);
581 void wallcycle_print(FILE *fplog, int nnodes, int npme, double realtime,
582 gmx_wallcycle_t wc, wallclock_gpu_t *gpu_t)
585 double c2t, tot, tot_gpu, tot_cpu_overlap, gpu_cpu_ratio, sum, tot_k;
586 int i, j, npp, nth_pp, nth_pme, nth_tot;
588 const char *hline = "-----------------------------------------------------------------------------";
595 nth_pp = wc->nthreads_pp;
596 nth_pme = wc->nthreads_pme;
598 cycles = wc->cycles_sum;
609 nth_tot = npp*nth_pp + npme*nth_pme;
611 tot = cycles[ewcRUN];
613 /* Conversion factor from cycles to seconds */
623 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");
625 fprintf(fplog, " Computing: Nodes Th. Count Wall t (s) G-Cycles %c\n", '%');
626 fprintf(fplog, "%s\n", hline);
628 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
630 if (!is_pme_subcounter(i))
632 print_cycles(fplog, c2t, wcn[i], nth_tot,
633 is_pme_counter(i) ? npme : npp,
634 is_pme_counter(i) ? nth_pme : nth_pp,
635 wc->wcc[i].n, cycles[i], tot);
639 if (wc->wcc_all != NULL)
641 for (i = 0; i < ewcNR; i++)
643 for (j = 0; j < ewcNR; j++)
645 snprintf(buf, 9, "%-9s", wcn[i]);
647 snprintf(buf+10, 9, "%-9s", wcn[j]);
649 print_cycles(fplog, c2t, buf, nth_tot,
650 is_pme_counter(i) ? npme : npp,
651 is_pme_counter(i) ? nth_pme : nth_pp,
652 wc->wcc_all[i*ewcNR+j].n,
653 wc->wcc_all[i*ewcNR+j].c,
658 print_cycles(fplog, c2t, "Rest", nth_tot, npp, -1, 0, tot-sum, tot);
659 fprintf(fplog, "%s\n", hline);
660 print_cycles(fplog, c2t, "Total", nth_tot, nnodes, -1, 0, tot, tot);
661 fprintf(fplog, "%s\n", hline);
663 if (wc->wcc[ewcPMEMESH].n > 0)
665 fprintf(fplog, "%s\n", hline);
666 for (i = ewcPPDURINGPME+1; i < ewcNR; i++)
668 if (is_pme_subcounter(i))
670 print_cycles(fplog, c2t, wcn[i], nth_tot,
671 is_pme_counter(i) ? npme : npp,
672 is_pme_counter(i) ? nth_pme : nth_pp,
673 wc->wcc[i].n, cycles[i], tot);
676 fprintf(fplog, "%s\n", hline);
679 #ifdef GMX_CYCLE_SUBCOUNTERS
680 fprintf(fplog, "%s\n", hline);
681 for (i = 0; i < ewcsNR; i++)
683 print_cycles(fplog, c2t, wcsn[i], nth_tot, npp, nth_pp,
684 wc->wcsc[i].n, cycles[ewcNR+i], tot);
686 fprintf(fplog, "%s\n", hline);
689 /* print GPU timing summary */
692 const char *k_log_str[2][2] = {
693 {"Nonbonded F kernel", "Nonbonded F+ene k."},
694 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
697 tot_gpu = gpu_t->pl_h2d_t + gpu_t->nb_h2d_t + gpu_t->nb_d2h_t;
699 /* add up the kernel timings */
701 for (i = 0; i < 2; i++)
703 for (j = 0; j < 2; j++)
705 tot_k += gpu_t->ktime[i][j].t;
710 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
711 if (wc->wcc[ewcPMEMESH].n > 0)
713 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
715 tot_cpu_overlap *= c2t * 1000; /* convert s to ms */
717 fprintf(fplog, "\n GPU timings\n%s\n", hline);
718 fprintf(fplog, " Computing: Count Wall t (s) ms/step %c\n", '%');
719 fprintf(fplog, "%s\n", hline);
720 print_gputimes(fplog, "Pair list H2D",
721 gpu_t->pl_h2d_c, gpu_t->pl_h2d_t, tot_gpu);
722 print_gputimes(fplog, "X / q H2D",
723 gpu_t->nb_c, gpu_t->nb_h2d_t, tot_gpu);
725 for (i = 0; i < 2; i++)
727 for (j = 0; j < 2; j++)
729 if (gpu_t->ktime[i][j].c)
731 print_gputimes(fplog, k_log_str[i][j],
732 gpu_t->ktime[i][j].c, gpu_t->ktime[i][j].t, tot_gpu);
737 print_gputimes(fplog, "F D2H", gpu_t->nb_c, gpu_t->nb_d2h_t, tot_gpu);
738 fprintf(fplog, "%s\n", hline);
739 print_gputimes(fplog, "Total ", gpu_t->nb_c, tot_gpu, tot_gpu);
740 fprintf(fplog, "%s\n", hline);
742 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
743 fprintf(fplog, "\nForce evaluation time GPU/CPU: %.3f ms/%.3f ms = %.3f\n",
744 tot_gpu/gpu_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
747 /* only print notes related to CPU-GPU load balance with PME */
748 if (wc->wcc[ewcPMEMESH].n > 0)
750 fprintf(fplog, "For optimal performance this ratio should be close to 1!\n");
752 /* print note if the imbalance is high with PME case in which
753 * CPU-GPU load balancing is possible */
754 if (gpu_cpu_ratio < 0.75 || gpu_cpu_ratio > 1.2)
756 /* Only the sim master calls this function, so always print to stderr */
757 if (gpu_cpu_ratio < 0.75)
761 /* The user could have used -notunepme,
762 * but we currently can't check that here.
764 md_print_warn(NULL, fplog,
765 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
766 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
767 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
771 /* We should not end up here, unless the box is
772 * too small for increasing the cut-off for PME tuning.
774 md_print_warn(NULL, fplog,
775 "\nNOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
776 " performance loss.");
779 if (gpu_cpu_ratio > 1.2)
781 md_print_warn(NULL, fplog,
782 "\nNOTE: The GPU has >20%% more load than the CPU. This imbalance causes\n"
783 " performance loss, consider using a shorter cut-off and a finer PME grid.");
789 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
790 (cycles[ewcDOMDEC] > tot*0.1 ||
791 cycles[ewcNS] > tot*0.1))
793 /* Only the sim master calls this function, so always print to stderr */
794 if (wc->wcc[ewcDOMDEC].n == 0)
796 md_print_warn(NULL, fplog,
797 "NOTE: %d %% of the run time was spent in pair search,\n"
798 " you might want to increase nstlist (this has no effect on accuracy)\n",
799 (int)(100*cycles[ewcNS]/tot+0.5));
803 md_print_warn(NULL, fplog,
804 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
805 " %d %% of the run time was spent in pair search,\n"
806 " you might want to increase nstlist (this has no effect on accuracy)\n",
807 (int)(100*cycles[ewcDOMDEC]/tot+0.5),
808 (int)(100*cycles[ewcNS]/tot+0.5));
812 if (cycles[ewcMoveE] > tot*0.05)
814 /* Only the sim master calls this function, so always print to stderr */
815 md_print_warn(NULL, fplog,
816 "NOTE: %d %% of the run time was spent communicating energies,\n"
817 " you might want to use the -gcom option of mdrun\n",
818 (int)(100*cycles[ewcMoveE]/tot+0.5));
822 extern gmx_large_int_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
829 return wc->reset_counters;
832 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_large_int_t reset_counters)
839 wc->reset_counters = reset_counters;
842 #ifdef GMX_CYCLE_SUBCOUNTERS
844 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
848 wc->wcsc[ewcs].start = gmx_cycles_read();
852 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
856 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
861 #endif /* GMX_CYCLE_SUBCOUNTERS */