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"
56 /* DEBUG_WCYCLE adds consistency checking for the counters.
57 * It checks if you stop a counter different from the last
58 * one that was opened and if you do nest too deep.
60 /* #define DEBUG_WCYCLE */
70 typedef struct gmx_wallcycle
73 /* variables for testing/debugging */
79 int counterlist[DEPTH_MAX];
83 gmx_cycles_t cycle_prev;
84 gmx_large_int_t reset_counters;
86 MPI_Comm mpi_comm_mygroup;
90 #ifdef GMX_CYCLE_SUBCOUNTERS
96 /* Each name should not exceed 19 characters */
97 static const char *wcn[ewcNR] =
98 { "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
99 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
100 "Comm. coord.", "Born radii", "Force", "Wait + Comm. F", "PME mesh",
101 "PME redist. X/F", "PME spread/gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve",
102 "PME wait for PP", "Wait + Recv. PME F", "Wait GPU nonlocal", "Wait GPU local", "NB X/F buffer ops.",
103 "Vsite spread", "Write traj.", "Update", "Constraints", "Comm. energies",
104 "Enforced rotation", "Add rot. forces", "Test" };
106 static const char *wcsn[ewcsNR] =
107 { "DD redist.", "DD NS grid + sort", "DD setup comm.",
108 "DD make top.", "DD make constr.", "DD top. other",
109 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
110 "Bonded F", "Nonbonded F", "Ewald F correction",
111 "NB X buffer ops.", "NB F buffer ops."
114 gmx_bool wallcycle_have_counter(void)
116 return gmx_cycles_have_counter();
119 gmx_wallcycle_t wallcycle_init(FILE *fplog,int resetstep,t_commrec *cr,
120 int nthreads_pp, int nthreads_pme)
125 if (!wallcycle_have_counter())
132 wc->wc_barrier = FALSE;
136 wc->reset_counters = resetstep;
137 wc->nthreads_pp = nthreads_pp;
138 wc->nthreads_pme = nthreads_pme;
139 wc->cycles_sum = NULL;
142 if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != NULL)
146 fprintf(fplog,"\nWill call MPI_Barrier before each cycle start/stop call\n\n");
148 wc->wc_barrier = TRUE;
149 wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
154 if (getenv("GMX_CYCLE_ALL") != NULL)
158 fprintf(fplog,"\nWill time all the code during the run\n\n");
160 snew(wc->wcc_all,ewcNR*ewcNR);
163 #ifdef GMX_CYCLE_SUBCOUNTERS
164 snew(wc->wcsc,ewcsNR);
174 void wallcycle_destroy(gmx_wallcycle_t wc)
185 if (wc->wcc_all != NULL)
189 #ifdef GMX_CYCLE_SUBCOUNTERS
190 if (wc->wcsc != NULL)
198 static void wallcycle_all_start(gmx_wallcycle_t wc,int ewc,gmx_cycles_t cycle)
201 wc->cycle_prev = cycle;
204 static void wallcycle_all_stop(gmx_wallcycle_t wc,int ewc,gmx_cycles_t cycle)
206 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
207 wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
212 static void debug_start_check(gmx_wallcycle_t wc, int ewc)
214 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
216 if (wc->count_depth < 0 || wc->count_depth >= DEPTH_MAX)
218 gmx_fatal(FARGS,"wallcycle counter depth out of range: %d",
221 wc->counterlist[wc->count_depth] = ewc;
225 static void debug_stop_check(gmx_wallcycle_t wc, int ewc)
229 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
231 if (wc->count_depth < 0)
233 gmx_fatal(FARGS,"wallcycle counter depth out of range when stopping %s: %d",wcn[ewc],wc->count_depth);
235 if (wc->counterlist[wc->count_depth] != ewc)
237 gmx_fatal(FARGS,"wallcycle mismatch at stop, start %s, stop %s",
238 wcn[wc->counterlist[wc->count_depth]],wcn[ewc]);
243 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
255 MPI_Barrier(wc->mpi_comm_mygroup);
260 debug_start_check(wc,ewc);
263 cycle = gmx_cycles_read();
264 wc->wcc[ewc].start = cycle;
265 if (wc->wcc_all != NULL)
270 wallcycle_all_start(wc,ewc,cycle);
272 else if (wc->wc_depth == 3)
274 wallcycle_all_stop(wc,ewc,cycle);
279 void wallcycle_start_nocount(gmx_wallcycle_t wc, int ewc)
286 wallcycle_start(wc, ewc);
290 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
292 gmx_cycles_t cycle,last;
302 MPI_Barrier(wc->mpi_comm_mygroup);
307 debug_stop_check(wc,ewc);
310 cycle = gmx_cycles_read();
311 last = cycle - wc->wcc[ewc].start;
312 wc->wcc[ewc].c += last;
319 wallcycle_all_stop(wc,ewc,cycle);
321 else if (wc->wc_depth == 2)
323 wallcycle_all_start(wc,ewc,cycle);
330 void wallcycle_reset_all(gmx_wallcycle_t wc)
339 for(i=0; i<ewcNR; i++)
346 for(i=0; i<ewcNR*ewcNR; i++)
348 wc->wcc_all[i].n = 0;
349 wc->wcc_all[i].c = 0;
352 #ifdef GMX_CYCLE_SUBCOUNTERS
353 for (i=0; i<ewcsNR; i++)
361 static gmx_bool is_pme_counter(int ewc)
363 return (ewc >= ewcPMEMESH && ewc <= ewcPMEWAITCOMM);
366 static gmx_bool is_pme_subcounter(int ewc)
368 return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
371 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
375 double cycles_n[ewcNR+ewcsNR],buf[ewcNR+ewcsNR],*cyc_all,*buf_all;
384 snew(wc->cycles_sum,ewcNR+ewcsNR);
385 cycles = wc->cycles_sum;
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,buf,nsum,MPI_DOUBLE,MPI_SUM,
483 for(i=0; i<nsum; i++)
488 if (wc->wcc_all != NULL)
490 snew(cyc_all,ewcNR*ewcNR);
491 snew(buf_all,ewcNR*ewcNR);
492 for(i=0; i<ewcNR*ewcNR; i++)
494 cyc_all[i] = wc->wcc_all[i].c;
496 MPI_Allreduce(cyc_all,buf_all,ewcNR*ewcNR,MPI_DOUBLE,MPI_SUM,
498 for(i=0; i<ewcNR*ewcNR; i++)
500 wc->wcc_all[i].c = buf_all[i];
509 static void print_cycles(FILE *fplog, double c2t, const char *name,
510 int nnodes_tot,int nnodes, int nthreads,
511 int n, double c, double tot)
521 snprintf(num,sizeof(num),"%10d",n);
523 snprintf(thstr, sizeof(thstr), "N/A");
525 snprintf(thstr, sizeof(thstr), "%4d", nthreads);
532 wallt = c*c2t*nnodes_tot/(double)nnodes;
533 fprintf(fplog," %-19s %4d %4s %10s %10.3f %12.3f %5.1f\n",
534 name,nnodes,thstr,num,wallt,c*1e-9,100*c/tot);
538 static void print_gputimes(FILE *fplog, const char *name,
539 int n, double t, double tot_t)
546 snprintf(num, sizeof(num), "%10d", n);
547 snprintf(avg_perf, sizeof(avg_perf), "%10.3f", t/n);
552 sprintf(avg_perf," ");
556 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
557 name, num, t/1000, avg_perf, 100 * t/tot_t);
561 fprintf(fplog, " %-29s %10s%12.3f %s %5.1f\n",
562 name, "", t/1000, avg_perf, 100.0);
566 void wallcycle_print(FILE *fplog, int nnodes, int npme, double realtime,
567 gmx_wallcycle_t wc, wallclock_gpu_t *gpu_t)
570 double c2t,tot,tot_gpu,tot_cpu_overlap,gpu_cpu_ratio,sum,tot_k;
571 int i,j,npp,nth_pp,nth_pme;
573 const char *hline = "-----------------------------------------------------------------------------";
580 nth_pp = wc->nthreads_pp;
581 nth_pme = wc->nthreads_pme;
583 cycles = wc->cycles_sum;
594 tot = cycles[ewcRUN];
596 /* Conversion factor from cycles to seconds */
606 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");
608 fprintf(fplog," Computing: Nodes Th. Count Wall t (s) G-Cycles %c\n",'%');
609 fprintf(fplog,"%s\n",hline);
611 for(i=ewcPPDURINGPME+1; i<ewcNR; i++)
613 if (!is_pme_subcounter(i))
615 print_cycles(fplog,c2t,wcn[i],nnodes,
616 is_pme_counter(i) ? npme : npp,
617 is_pme_counter(i) ? nth_pme : nth_pp,
618 wc->wcc[i].n,cycles[i],tot);
622 if (wc->wcc_all != NULL)
624 for(i=0; i<ewcNR; i++)
626 for(j=0; j<ewcNR; j++)
628 snprintf(buf,9,"%-9s",wcn[i]);
630 snprintf(buf+10,9,"%-9s",wcn[j]);
632 print_cycles(fplog,c2t,buf,nnodes,
633 is_pme_counter(i) ? npme : npp,
634 is_pme_counter(i) ? nth_pme : nth_pp,
635 wc->wcc_all[i*ewcNR+j].n,
636 wc->wcc_all[i*ewcNR+j].c,
641 print_cycles(fplog,c2t,"Rest",npp,npp,-1,0,tot-sum,tot);
642 fprintf(fplog,"%s\n",hline);
643 print_cycles(fplog,c2t,"Total",nnodes,nnodes,-1,0,tot,tot);
644 fprintf(fplog,"%s\n",hline);
646 if (wc->wcc[ewcPMEMESH].n > 0)
648 fprintf(fplog,"%s\n",hline);
649 for(i=ewcPPDURINGPME+1; i<ewcNR; i++)
651 if (is_pme_subcounter(i))
653 print_cycles(fplog,c2t,wcn[i],nnodes,
654 is_pme_counter(i) ? npme : npp,
655 is_pme_counter(i) ? nth_pme : nth_pp,
656 wc->wcc[i].n,cycles[i],tot);
659 fprintf(fplog,"%s\n",hline);
662 #ifdef GMX_CYCLE_SUBCOUNTERS
663 fprintf(fplog,"%s\n",hline);
664 for(i=0; i<ewcsNR; i++)
666 print_cycles(fplog,c2t,wcsn[i],nnodes,npp,nth_pp,
667 wc->wcsc[i].n,cycles[ewcNR+i],tot);
669 fprintf(fplog,"%s\n",hline);
672 /* print GPU timing summary */
675 const char *k_log_str[2][2] = {
676 {"Nonbonded F kernel", "Nonbonded F+ene k."},
677 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}};
679 tot_gpu = gpu_t->pl_h2d_t + gpu_t->nb_h2d_t + gpu_t->nb_d2h_t;
681 /* add up the kernel timings */
683 for (i = 0; i < 2; i++)
685 for(j = 0; j < 2; j++)
687 tot_k += gpu_t->ktime[i][j].t;
692 tot_cpu_overlap = wc->wcc[ewcFORCE].c;
693 if (wc->wcc[ewcPMEMESH].n > 0)
695 tot_cpu_overlap += wc->wcc[ewcPMEMESH].c;
697 tot_cpu_overlap *= c2t * 1000; /* convert s to ms */
699 fprintf(fplog, "\n GPU timings\n%s\n", hline);
700 fprintf(fplog," Computing: Count Wall t (s) ms/step %c\n",'%');
701 fprintf(fplog, "%s\n", hline);
702 print_gputimes(fplog, "Pair list H2D",
703 gpu_t->pl_h2d_c, gpu_t->pl_h2d_t, tot_gpu);
704 print_gputimes(fplog, "X / q H2D",
705 gpu_t->nb_c, gpu_t->nb_h2d_t, tot_gpu);
707 for (i = 0; i < 2; i++)
709 for(j = 0; j < 2; j++)
711 if (gpu_t->ktime[i][j].c)
713 print_gputimes(fplog, k_log_str[i][j],
714 gpu_t->ktime[i][j].c, gpu_t->ktime[i][j].t, tot_gpu);
719 print_gputimes(fplog, "F D2H", gpu_t->nb_c, gpu_t->nb_d2h_t, tot_gpu);
720 fprintf(fplog, "%s\n", hline);
721 print_gputimes(fplog, "Total ", gpu_t->nb_c, tot_gpu, tot_gpu);
722 fprintf(fplog, "%s\n", hline);
724 gpu_cpu_ratio = tot_gpu/tot_cpu_overlap;
725 fprintf(fplog, "\n Force evaluation time GPU/CPU: %.3f ms/%.3f ms = %.3f\n",
726 tot_gpu/gpu_t->nb_c, tot_cpu_overlap/wc->wcc[ewcFORCE].n,
729 /* only print notes related to CPU-GPU load balance with PME */
730 if (wc->wcc[ewcPMEMESH].n > 0)
732 fprintf(fplog, "For optimal performance this ratio should be close to 1!\n");
734 /* print note if the imbalance is high with PME case in which
735 * CPU-GPU load balancing is possible */
736 if (gpu_cpu_ratio < 0.75 || gpu_cpu_ratio > 1.2)
738 /* Only the sim master calls this function, so always print to stderr */
739 if (gpu_cpu_ratio < 0.75)
743 /* The user could have used -notunepme,
744 * but we currently can't check that here.
746 md_print_warn(NULL,fplog,
747 "NOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
748 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
749 " In that case, try setting the DD grid manually (-dd) or lowering -dds.\n");
753 /* We should not end up here, unless the box is
754 * too small for increasing the cut-off for PME tuning.
756 md_print_warn(NULL,fplog,
757 "NOTE: The GPU has >25%% less load than the CPU. This imbalance causes\n"
758 " performance loss.\n");
761 if (gpu_cpu_ratio > 1.2)
763 md_print_warn(NULL,fplog,
764 "NOTE: The GPU has >20%% more load than the CPU. This imbalance causes\n"
765 " performance loss, consider using a shorter cut-off and a finer PME grid.\n");
771 if (wc->wcc[ewcNB_XF_BUF_OPS].n > 0 &&
772 (cycles[ewcDOMDEC] > tot*0.1 ||
773 cycles[ewcNS] > tot*0.1))
775 /* Only the sim master calls this function, so always print to stderr */
776 if (wc->wcc[ewcDOMDEC].n == 0)
778 md_print_warn(NULL,fplog,
779 "NOTE: %d %% of the run time was spent in pair search,\n"
780 " you might want to increase nstlist (this has no effect on accuracy)\n",
781 (int)(100*cycles[ewcNS]/tot+0.5));
785 md_print_warn(NULL,fplog,
786 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
787 " %d %% of the run time was spent in pair search,\n"
788 " you might want to increase nstlist (this has no effect on accuracy)\n",
789 (int)(100*cycles[ewcDOMDEC]/tot+0.5),
790 (int)(100*cycles[ewcNS]/tot+0.5));
794 if (cycles[ewcMoveE] > tot*0.05)
796 /* Only the sim master calls this function, so always print to stderr */
797 md_print_warn(NULL,fplog,
798 "NOTE: %d %% of the run time was spent communicating energies,\n"
799 " you might want to use the -gcom option of mdrun\n",
800 (int)(100*cycles[ewcMoveE]/tot+0.5));
804 extern gmx_large_int_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
811 return wc->reset_counters;
814 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_large_int_t reset_counters)
819 wc->reset_counters = reset_counters;
822 #ifdef GMX_CYCLE_SUBCOUNTERS
824 void wallcycle_sub_start(gmx_wallcycle_t wc, int ewcs)
828 wc->wcsc[ewcs].start = gmx_cycles_read();
832 void wallcycle_sub_stop(gmx_wallcycle_t wc, int ewcs)
836 wc->wcsc[ewcs].c += gmx_cycles_read() - wc->wcsc[ewcs].start;
841 #endif /* GMX_CYCLE_SUBCOUNTERS */