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,2015,2016,2017, 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 "wallcycle.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/timing/cyclecounter.h"
49 #include "gromacs/timing/gpu_timing.h"
50 #include "gromacs/timing/wallcyclereporting.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/gmxassert.h"
53 #include "gromacs/utility/gmxmpi.h"
54 #include "gromacs/utility/logger.h"
55 #include "gromacs/utility/smalloc.h"
56 #include "gromacs/utility/snprintf.h"
58 static const bool useCycleSubcounters
= GMX_CYCLE_SUBCOUNTERS
;
60 /* DEBUG_WCYCLE adds consistency checking for the counters.
61 * It checks if you stop a counter different from the last
62 * one that was opened and if you do nest too deep.
64 /* #define DEBUG_WCYCLE */
67 #include "gromacs/utility/fatalerror.h"
77 typedef struct gmx_wallcycle
80 /* did we detect one or more invalid cycle counts */
81 gmx_bool haveInvalidCount
;
82 /* variables for testing/debugging */
88 int counterlist
[DEPTH_MAX
];
92 gmx_cycles_t cycle_prev
;
93 gmx_int64_t reset_counters
;
95 MPI_Comm mpi_comm_mygroup
;
100 /* Each name should not exceed 19 printing characters
101 (ie. terminating null can be twentieth) */
102 static const char *wcn
[ewcNR
] =
104 "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load",
105 "DD comm. bounds", "Vsite constr.", "Send X to PME", "Neighbor search", "Launch GPU ops.",
106 "Comm. coord.", "Born radii", "Force", "Wait + Comm. F", "PME mesh",
107 "PME redist. X/F", "PME spread", "PME gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve LJ", "PME solve Elec",
108 "PME wait for PP", "Wait + Recv. PME F",
109 "Wait PME GPU spread", "Wait PME GPU gather", "Reduce GPU PME F",
110 "Wait GPU NB nonloc.", "Wait GPU NB local", "NB X/F buffer ops.",
111 "Vsite spread", "COM pull force",
112 "Write traj.", "Update", "Constraints", "Comm. energies",
113 "Enforced rotation", "Add rot. forces", "Position swapping", "IMD", "Test"
116 static const char *wcsn
[ewcsNR
] =
118 "DD redist.", "DD NS grid + sort", "DD setup comm.",
119 "DD make top.", "DD make constr.", "DD top. other",
120 "NS grid local", "NS grid non-loc.", "NS search local", "NS search non-loc.",
124 "Listed buffer ops.",
127 "Launch NB GPU tasks",
128 "Launch PME GPU tasks",
129 "Ewald F correction",
134 /* PME GPU timing events' names - correspond to the enum in the gpu_timing.h */
135 static const char *PMEStageNames
[] =
146 gmx_bool
wallcycle_have_counter(void)
148 return gmx_cycles_have_counter();
151 gmx_wallcycle_t
wallcycle_init(FILE *fplog
, int resetstep
, t_commrec gmx_unused
*cr
)
156 if (!wallcycle_have_counter())
163 wc
->haveInvalidCount
= FALSE
;
164 wc
->wc_barrier
= FALSE
;
165 wc
->wcc_all
= nullptr;
168 wc
->reset_counters
= resetstep
;
171 if (PAR(cr
) && getenv("GMX_CYCLE_BARRIER") != nullptr)
175 fprintf(fplog
, "\nWill call MPI_Barrier before each cycle start/stop call\n\n");
177 wc
->wc_barrier
= TRUE
;
178 wc
->mpi_comm_mygroup
= cr
->mpi_comm_mygroup
;
182 snew(wc
->wcc
, ewcNR
);
183 if (getenv("GMX_CYCLE_ALL") != nullptr)
187 fprintf(fplog
, "\nWill time all the code during the run\n\n");
189 snew(wc
->wcc_all
, ewcNR
*ewcNR
);
192 if (useCycleSubcounters
)
194 snew(wc
->wcsc
, ewcsNR
);
204 /* TODO: Should be called from finish_run() or runner()
205 void wallcycle_destroy(gmx_wallcycle_t wc)
212 if (wc->wcc != nullptr)
216 if (wc->wcc_all != nullptr)
220 if (wc->wcsc != nullptr)
228 static void wallcycle_all_start(gmx_wallcycle_t wc
, int ewc
, gmx_cycles_t cycle
)
231 wc
->cycle_prev
= cycle
;
234 static void wallcycle_all_stop(gmx_wallcycle_t wc
, int ewc
, gmx_cycles_t cycle
)
236 wc
->wcc_all
[wc
->ewc_prev
*ewcNR
+ewc
].n
+= 1;
237 wc
->wcc_all
[wc
->ewc_prev
*ewcNR
+ewc
].c
+= cycle
- wc
->cycle_prev
;
242 static void debug_start_check(gmx_wallcycle_t wc
, int ewc
)
244 /* fprintf(stderr,"wcycle_start depth %d, %s\n",wc->count_depth,wcn[ewc]); */
246 if (wc
->count_depth
< 0 || wc
->count_depth
>= DEPTH_MAX
)
248 gmx_fatal(FARGS
, "wallcycle counter depth out of range: %d",
251 wc
->counterlist
[wc
->count_depth
] = ewc
;
255 static void debug_stop_check(gmx_wallcycle_t wc
, int ewc
)
259 /* fprintf(stderr,"wcycle_stop depth %d, %s\n",wc->count_depth,wcn[ewc]); */
261 if (wc
->count_depth
< 0)
263 gmx_fatal(FARGS
, "wallcycle counter depth out of range when stopping %s: %d", wcn
[ewc
], wc
->count_depth
);
265 if (wc
->counterlist
[wc
->count_depth
] != ewc
)
267 gmx_fatal(FARGS
, "wallcycle mismatch at stop, start %s, stop %s",
268 wcn
[wc
->counterlist
[wc
->count_depth
]], wcn
[ewc
]);
273 void wallcycle_start(gmx_wallcycle_t wc
, int ewc
)
285 MPI_Barrier(wc
->mpi_comm_mygroup
);
290 debug_start_check(wc
, ewc
);
293 cycle
= gmx_cycles_read();
294 wc
->wcc
[ewc
].start
= cycle
;
295 if (wc
->wcc_all
!= nullptr)
300 wallcycle_all_start(wc
, ewc
, cycle
);
302 else if (wc
->wc_depth
== 3)
304 wallcycle_all_stop(wc
, ewc
, cycle
);
309 void wallcycle_start_nocount(gmx_wallcycle_t wc
, int ewc
)
316 wallcycle_start(wc
, ewc
);
320 double wallcycle_stop(gmx_wallcycle_t wc
, int ewc
)
322 gmx_cycles_t cycle
, last
;
332 MPI_Barrier(wc
->mpi_comm_mygroup
);
337 debug_stop_check(wc
, ewc
);
340 /* When processes or threads migrate between cores, the cycle counting
341 * can get messed up if the cycle counter on different cores are not
342 * synchronized. When this happens we expect both large negative and
343 * positive cycle differences. We can detect negative cycle differences.
344 * Detecting too large positive counts if difficult, since count can be
345 * large, especially for ewcRUN. If we detect a negative count,
346 * we will not print the cycle accounting table.
348 cycle
= gmx_cycles_read();
349 if (cycle
>= wc
->wcc
[ewc
].start
)
351 last
= cycle
- wc
->wcc
[ewc
].start
;
356 wc
->haveInvalidCount
= TRUE
;
358 wc
->wcc
[ewc
].c
+= last
;
365 wallcycle_all_stop(wc
, ewc
, cycle
);
367 else if (wc
->wc_depth
== 2)
369 wallcycle_all_start(wc
, ewc
, cycle
);
376 void wallcycle_get(gmx_wallcycle_t wc
, int ewc
, int *n
, double *c
)
379 *c
= static_cast<double>(wc
->wcc
[ewc
].c
);
382 void wallcycle_reset_all(gmx_wallcycle_t wc
)
391 for (i
= 0; i
< ewcNR
; i
++)
396 wc
->haveInvalidCount
= FALSE
;
400 for (i
= 0; i
< ewcNR
*ewcNR
; i
++)
402 wc
->wcc_all
[i
].n
= 0;
403 wc
->wcc_all
[i
].c
= 0;
408 for (i
= 0; i
< ewcsNR
; i
++)
416 static gmx_bool
is_pme_counter(int ewc
)
418 return (ewc
>= ewcPMEMESH
&& ewc
<= ewcPMEWAITCOMM
);
421 static gmx_bool
is_pme_subcounter(int ewc
)
423 return (ewc
>= ewcPME_REDISTXF
&& ewc
< ewcPMEWAITCOMM
);
426 /* Subtract counter ewc_sub timed inside a timing block for ewc_main */
427 static void subtract_cycles(wallcc_t
*wcc
, int ewc_main
, int ewc_sub
)
429 if (wcc
[ewc_sub
].n
> 0)
431 if (wcc
[ewc_main
].c
>= wcc
[ewc_sub
].c
)
433 wcc
[ewc_main
].c
-= wcc
[ewc_sub
].c
;
437 /* Something is wrong with the cycle counting */
443 void wallcycle_scale_by_num_threads(gmx_wallcycle_t wc
, bool isPmeRank
, int nthreads_pp
, int nthreads_pme
)
450 for (int i
= 0; i
< ewcNR
; i
++)
452 if (is_pme_counter(i
) || (i
== ewcRUN
&& isPmeRank
))
454 wc
->wcc
[i
].c
*= nthreads_pme
;
458 for (int j
= 0; j
< ewcNR
; j
++)
460 wc
->wcc_all
[i
*ewcNR
+j
].c
*= nthreads_pme
;
466 wc
->wcc
[i
].c
*= nthreads_pp
;
470 for (int j
= 0; j
< ewcNR
; j
++)
472 wc
->wcc_all
[i
*ewcNR
+j
].c
*= nthreads_pp
;
477 if (useCycleSubcounters
&& wc
->wcsc
&& !isPmeRank
)
479 for (int i
= 0; i
< ewcsNR
; i
++)
481 wc
->wcsc
[i
].c
*= nthreads_pp
;
486 /* TODO Make an object for this function to return, containing some
487 * vectors of something like wallcc_t for the summed wcc, wcc_all and
488 * wcsc, AND the original wcc for rank 0.
490 * The GPU timing is reported only for rank 0, so we want to preserve
491 * the original wcycle on that rank. Rank 0 also reports the global
492 * counts before that, so needs something to contain the global data
493 * without over-writing the rank-0 data. The current implementation
494 * uses cycles_sum to manage this, which works OK now because wcsc and
495 * wcc_all are unused by the GPU reporting, but it is not satisfactory
496 * for the future. Also, there's no need for MPI_Allreduce, since
497 * only MASTERRANK uses any of the results. */
498 WallcycleCounts
wallcycle_sum(t_commrec
*cr
, gmx_wallcycle_t wc
)
500 WallcycleCounts cycles_sum
;
502 double cycles
[ewcNR
+ewcsNR
];
504 double cycles_n
[ewcNR
+ewcsNR
+1];
511 /* Default construction of std::array of non-class T can leave
512 the values indeterminate, just like a C array, and icc
520 subtract_cycles(wcc
, ewcDOMDEC
, ewcDDCOMMLOAD
);
521 subtract_cycles(wcc
, ewcDOMDEC
, ewcDDCOMMBOUND
);
523 subtract_cycles(wcc
, ewcPME_FFT
, ewcPME_FFTCOMM
);
525 if (cr
->npmenodes
== 0)
527 /* All nodes do PME (or no PME at all) */
528 subtract_cycles(wcc
, ewcFORCE
, ewcPMEMESH
);
532 /* The are PME-only nodes */
533 if (wcc
[ewcPMEMESH
].n
> 0)
535 /* This must be a PME only node, calculate the Wait + Comm. time */
536 GMX_ASSERT(wcc
[ewcRUN
].c
>= wcc
[ewcPMEMESH
].c
, "Total run ticks must be greater than PME-only ticks");
537 wcc
[ewcPMEWAITCOMM
].c
= wcc
[ewcRUN
].c
- wcc
[ewcPMEMESH
].c
;
541 /* Store the cycles in a double buffer for summing */
542 for (i
= 0; i
< ewcNR
; i
++)
545 cycles_n
[i
] = static_cast<double>(wcc
[i
].n
);
547 cycles
[i
] = static_cast<double>(wcc
[i
].c
);
552 for (i
= 0; i
< ewcsNR
; i
++)
555 cycles_n
[ewcNR
+i
] = static_cast<double>(wc
->wcsc
[i
].n
);
557 cycles
[ewcNR
+i
] = static_cast<double>(wc
->wcsc
[i
].c
);
565 double buf
[ewcNR
+ewcsNR
+1];
567 // TODO this code is used only at the end of the run, so we
568 // can just do a simple reduce of haveInvalidCount in
569 // wallcycle_print, and avoid bugs
570 cycles_n
[nsum
] = (wc
->haveInvalidCount
> 0 ? 1 : 0);
571 // TODO Use MPI_Reduce
572 MPI_Allreduce(cycles_n
, buf
, nsum
+ 1, MPI_DOUBLE
, MPI_MAX
,
574 for (i
= 0; i
< ewcNR
; i
++)
576 wcc
[i
].n
= static_cast<int>(buf
[i
] + 0.5);
578 wc
->haveInvalidCount
= (buf
[nsum
] > 0);
581 for (i
= 0; i
< ewcsNR
; i
++)
583 wc
->wcsc
[i
].n
= static_cast<int>(buf
[ewcNR
+i
] + 0.5);
587 // TODO Use MPI_Reduce
588 MPI_Allreduce(cycles
, cycles_sum
.data(), nsum
, MPI_DOUBLE
, MPI_SUM
,
591 if (wc
->wcc_all
!= nullptr)
593 double *buf_all
, *cyc_all
;
595 snew(cyc_all
, ewcNR
*ewcNR
);
596 snew(buf_all
, ewcNR
*ewcNR
);
597 for (i
= 0; i
< ewcNR
*ewcNR
; i
++)
599 cyc_all
[i
] = wc
->wcc_all
[i
].c
;
601 // TODO Use MPI_Reduce
602 MPI_Allreduce(cyc_all
, buf_all
, ewcNR
*ewcNR
, MPI_DOUBLE
, MPI_SUM
,
604 for (i
= 0; i
< ewcNR
*ewcNR
; i
++)
606 wc
->wcc_all
[i
].c
= static_cast<gmx_cycles_t
>(buf_all
[i
]);
615 for (i
= 0; i
< nsum
; i
++)
617 cycles_sum
[i
] = cycles
[i
];
624 static void print_cycles(FILE *fplog
, double c2t
, const char *name
,
625 int nnodes
, int nthreads
,
626 int ncalls
, double c_sum
, double tot
)
628 char nnodes_str
[STRLEN
];
629 char nthreads_str
[STRLEN
];
630 char ncalls_str
[STRLEN
];
632 double percentage
= (tot
> 0.) ? (100. * c_sum
/ tot
) : 0.;
638 snprintf(ncalls_str
, sizeof(ncalls_str
), "%10d", ncalls
);
641 snprintf(nnodes_str
, sizeof(nnodes_str
), "N/A");
645 snprintf(nnodes_str
, sizeof(nnodes_str
), "%4d", nnodes
);
649 snprintf(nthreads_str
, sizeof(nthreads_str
), "N/A");
653 snprintf(nthreads_str
, sizeof(nthreads_str
), "%4d", nthreads
);
662 /* Convert the cycle count to wallclock time for this task */
665 fprintf(fplog
, " %-19.19s %4s %4s %10s %10.3f %14.3f %5.1f\n",
666 name
, nnodes_str
, nthreads_str
, ncalls_str
, wallt
,
667 c_sum
*1e-9, percentage
);
671 static void print_gputimes(FILE *fplog
, const char *name
,
672 int n
, double t
, double tot_t
)
679 snprintf(num
, sizeof(num
), "%10d", n
);
680 snprintf(avg_perf
, sizeof(avg_perf
), "%10.3f", t
/n
);
685 sprintf(avg_perf
, " ");
687 if (t
!= tot_t
&& tot_t
> 0)
689 fprintf(fplog
, " %-29s %10s%12.3f %s %5.1f\n",
690 name
, num
, t
/1000, avg_perf
, 100 * t
/tot_t
);
694 fprintf(fplog
, " %-29s %10s%12.3f %s %5.1f\n",
695 name
, "", t
/1000, avg_perf
, 100.0);
699 static void print_header(FILE *fplog
, int nrank_pp
, int nth_pp
, int nrank_pme
, int nth_pme
)
701 int nrank_tot
= nrank_pp
+ nrank_pme
;
704 fprintf(fplog
, "On %d MPI rank%s", nrank_tot
, nrank_tot
== 1 ? "" : "s");
707 fprintf(fplog
, ", each using %d OpenMP threads", nth_pp
);
709 /* Don't report doing PP+PME, because we can't tell here if
710 * this is RF, etc. */
714 fprintf(fplog
, "On %d MPI rank%s doing PP", nrank_pp
, nrank_pp
== 1 ? "" : "s");
717 fprintf(fplog
, ",%s using %d OpenMP threads", nrank_pp
> 1 ? " each" : "", nth_pp
);
719 fprintf(fplog
, ", and\non %d MPI rank%s doing PME", nrank_pme
, nrank_pme
== 1 ? "" : "s");
722 fprintf(fplog
, ",%s using %d OpenMP threads", nrank_pme
> 1 ? " each" : "", nth_pme
);
726 fprintf(fplog
, "\n\n");
727 fprintf(fplog
, " Computing: Num Num Call Wall time Giga-Cycles\n");
728 fprintf(fplog
, " Ranks Threads Count (s) total sum %%\n");
732 void wallcycle_print(FILE *fplog
, const gmx::MDLogger
&mdlog
, int nnodes
, int npme
,
733 int nth_pp
, int nth_pme
, double realtime
,
734 gmx_wallcycle_t wc
, const WallcycleCounts
&cyc_sum
,
735 const gmx_wallclock_gpu_nbnxn_t
*gpu_nbnxn_t
,
736 const gmx_wallclock_gpu_pme_t
*gpu_pme_t
)
738 double tot
, tot_for_pp
, tot_for_rest
, tot_cpu_overlap
, gpu_cpu_ratio
;
739 double c2t
, c2t_pp
, c2t_pme
= 0;
740 int i
, j
, npp
, nth_tot
;
742 const char *hline
= "-----------------------------------------------------------------------------";
749 GMX_ASSERT(nth_pp
> 0, "Number of particle-particle threads must be >0");
750 GMX_ASSERT(nth_pme
> 0, "Number of PME threads must be >0");
751 GMX_ASSERT(nnodes
> 0, "Number of nodes must be >0");
752 GMX_ASSERT(npme
>= 0, "Number of PME nodes cannot be negative");
754 /* npme is the number of PME-only ranks used, and we always do PP work */
755 GMX_ASSERT(npp
> 0, "Number of particle-particle nodes must be >0");
757 nth_tot
= npp
*nth_pp
+ npme
*nth_pme
;
759 /* When using PME-only nodes, the next line is valid for both
760 PP-only and PME-only nodes because they started ewcRUN at the
762 tot
= cyc_sum
[ewcRUN
];
767 /* TODO This is heavy handed, but until someone reworks the
768 code so that it is provably robust with respect to
769 non-positive values for all possible timer and cycle
770 counters, there is less value gained from printing whatever
771 timing data might still be sensible for some non-Jenkins
772 run, than is lost from diagnosing Jenkins FP exceptions on
773 runs about whose execution time we don't care. */
774 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
775 "WARNING: A total of %f CPU cycles was recorded, so mdrun cannot print a time accounting",
780 if (wc
->haveInvalidCount
)
782 GMX_LOG(mdlog
.warning
).asParagraph().appendText("NOTE: Detected invalid cycle counts, probably because threads moved between CPU cores that do not have synchronized cycle counters. Will not print the cycle accounting.");
787 /* Conversion factor from cycles to seconds */
789 c2t_pp
= c2t
* nth_tot
/ static_cast<double>(npp
*nth_pp
);
792 c2t_pme
= c2t
* nth_tot
/ static_cast<double>(npme
*nth_pme
);
799 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");
801 print_header(fplog
, npp
, nth_pp
, npme
, nth_pme
);
803 fprintf(fplog
, "%s\n", hline
);
804 for (i
= ewcPPDURINGPME
+1; i
< ewcNR
; i
++)
806 if (is_pme_subcounter(i
))
808 /* Do not count these at all */
810 else if (npme
> 0 && is_pme_counter(i
))
812 /* Print timing information for PME-only nodes, but add an
813 * asterisk so the reader of the table can know that the
814 * walltimes are not meant to add up. The asterisk still
815 * fits in the required maximum of 19 characters. */
817 snprintf(buffer
, STRLEN
, "%s *", wcn
[i
]);
818 print_cycles(fplog
, c2t_pme
, buffer
,
820 wc
->wcc
[i
].n
, cyc_sum
[i
], tot
);
824 /* Print timing information when it is for a PP or PP+PME
826 print_cycles(fplog
, c2t_pp
, wcn
[i
],
828 wc
->wcc
[i
].n
, cyc_sum
[i
], tot
);
829 tot_for_pp
+= cyc_sum
[i
];
832 if (wc
->wcc_all
!= nullptr)
834 for (i
= 0; i
< ewcNR
; i
++)
836 for (j
= 0; j
< ewcNR
; j
++)
838 snprintf(buf
, 20, "%-9.9s %-9.9s", wcn
[i
], wcn
[j
]);
839 print_cycles(fplog
, c2t_pp
, buf
,
841 wc
->wcc_all
[i
*ewcNR
+j
].n
,
842 wc
->wcc_all
[i
*ewcNR
+j
].c
,
847 tot_for_rest
= tot
* npp
* nth_pp
/ static_cast<double>(nth_tot
);
848 print_cycles(fplog
, c2t_pp
, "Rest",
850 -1, tot_for_rest
- tot_for_pp
, tot
);
851 fprintf(fplog
, "%s\n", hline
);
852 print_cycles(fplog
, c2t
, "Total",
855 fprintf(fplog
, "%s\n", hline
);
860 "(*) Note that with separate PME ranks, the walltime column actually sums to\n"
861 " twice the total reported, but the cycle count total and %% are correct.\n"
865 if (wc
->wcc
[ewcPMEMESH
].n
> 0)
867 fprintf(fplog
, " Breakdown of PME mesh computation\n");
868 fprintf(fplog
, "%s\n", hline
);
869 for (i
= ewcPPDURINGPME
+1; i
< ewcNR
; i
++)
871 if (is_pme_subcounter(i
))
873 print_cycles(fplog
, npme
> 0 ? c2t_pme
: c2t_pp
, wcn
[i
],
874 npme
> 0 ? npme
: npp
, nth_pme
,
875 wc
->wcc
[i
].n
, cyc_sum
[i
], tot
);
878 fprintf(fplog
, "%s\n", hline
);
881 if (useCycleSubcounters
&& wc
->wcsc
)
883 fprintf(fplog
, " Breakdown of PP computation\n");
884 fprintf(fplog
, "%s\n", hline
);
885 for (i
= 0; i
< ewcsNR
; i
++)
887 print_cycles(fplog
, c2t_pp
, wcsn
[i
],
889 wc
->wcsc
[i
].n
, cyc_sum
[ewcNR
+i
], tot
);
891 fprintf(fplog
, "%s\n", hline
);
894 /* print GPU timing summary */
895 double tot_gpu
= 0.0;
898 for (size_t k
= 0; k
< gtPME_EVENT_COUNT
; k
++)
900 tot_gpu
+= gpu_pme_t
->timing
[k
].t
;
905 const char *k_log_str
[2][2] = {
906 {"Nonbonded F kernel", "Nonbonded F+ene k."},
907 {"Nonbonded F+prune k.", "Nonbonded F+ene+prune k."}
909 tot_gpu
+= gpu_nbnxn_t
->pl_h2d_t
+ gpu_nbnxn_t
->nb_h2d_t
+ gpu_nbnxn_t
->nb_d2h_t
;
911 /* add up the kernel timings */
912 for (i
= 0; i
< 2; i
++)
914 for (j
= 0; j
< 2; j
++)
916 tot_gpu
+= gpu_nbnxn_t
->ktime
[i
][j
].t
;
919 tot_gpu
+= gpu_nbnxn_t
->pruneTime
.t
;
921 tot_cpu_overlap
= wc
->wcc
[ewcFORCE
].c
;
922 if (wc
->wcc
[ewcPMEMESH
].n
> 0)
924 tot_cpu_overlap
+= wc
->wcc
[ewcPMEMESH
].c
;
926 tot_cpu_overlap
*= realtime
*1000/tot
; /* convert s to ms */
928 fprintf(fplog
, "\n GPU timings\n%s\n", hline
);
929 fprintf(fplog
, " Computing: Count Wall t (s) ms/step %c\n", '%');
930 fprintf(fplog
, "%s\n", hline
);
931 print_gputimes(fplog
, "Pair list H2D",
932 gpu_nbnxn_t
->pl_h2d_c
, gpu_nbnxn_t
->pl_h2d_t
, tot_gpu
);
933 print_gputimes(fplog
, "X / q H2D",
934 gpu_nbnxn_t
->nb_c
, gpu_nbnxn_t
->nb_h2d_t
, tot_gpu
);
936 for (i
= 0; i
< 2; i
++)
938 for (j
= 0; j
< 2; j
++)
940 if (gpu_nbnxn_t
->ktime
[i
][j
].c
)
942 print_gputimes(fplog
, k_log_str
[i
][j
],
943 gpu_nbnxn_t
->ktime
[i
][j
].c
, gpu_nbnxn_t
->ktime
[i
][j
].t
, tot_gpu
);
949 for (size_t k
= 0; k
< gtPME_EVENT_COUNT
; k
++)
951 if (gpu_pme_t
->timing
[k
].c
)
953 print_gputimes(fplog
, PMEStageNames
[k
],
954 gpu_pme_t
->timing
[k
].c
,
955 gpu_pme_t
->timing
[k
].t
,
960 if (gpu_nbnxn_t
->pruneTime
.c
)
962 print_gputimes(fplog
, "Pruning kernel", gpu_nbnxn_t
->pruneTime
.c
, gpu_nbnxn_t
->pruneTime
.t
, tot_gpu
);
964 print_gputimes(fplog
, "F D2H", gpu_nbnxn_t
->nb_c
, gpu_nbnxn_t
->nb_d2h_t
, tot_gpu
);
965 fprintf(fplog
, "%s\n", hline
);
966 print_gputimes(fplog
, "Total ", gpu_nbnxn_t
->nb_c
, tot_gpu
, tot_gpu
);
967 fprintf(fplog
, "%s\n", hline
);
968 if (gpu_nbnxn_t
->dynamicPruneTime
.c
)
970 /* We print the dynamic pruning kernel timings after a separator
971 * and avoid adding it to tot_gpu as this is not in the force
972 * overlap. We print the fraction as relative to the rest.
974 print_gputimes(fplog
, "*Dynamic pruning", gpu_nbnxn_t
->dynamicPruneTime
.c
, gpu_nbnxn_t
->dynamicPruneTime
.t
, tot_gpu
);
975 fprintf(fplog
, "%s\n", hline
);
977 gpu_cpu_ratio
= tot_gpu
/tot_cpu_overlap
;
978 if (gpu_nbnxn_t
->nb_c
> 0 && wc
->wcc
[ewcFORCE
].n
> 0)
980 // FIXME the code below is not updated for PME on GPU
981 fprintf(fplog
, "\nAverage per-step force GPU/CPU evaluation time ratio: %.3f ms/%.3f ms = %.3f\n",
982 tot_gpu
/gpu_nbnxn_t
->nb_c
, tot_cpu_overlap
/wc
->wcc
[ewcFORCE
].n
,
986 /* only print notes related to CPU-GPU load balance with PME */
987 if (wc
->wcc
[ewcPMEMESH
].n
> 0)
989 fprintf(fplog
, "For optimal performance this ratio should be close to 1!\n");
991 /* print note if the imbalance is high with PME case in which
992 * CPU-GPU load balancing is possible */
993 if (gpu_cpu_ratio
< 0.75 || gpu_cpu_ratio
> 1.2)
995 /* Only the sim master calls this function, so always print to stderr */
996 if (gpu_cpu_ratio
< 0.75)
1000 /* The user could have used -notunepme,
1001 * but we currently can't check that here.
1003 GMX_LOG(mdlog
.warning
).asParagraph().appendText(
1004 "NOTE: The GPU has >25% less load than the CPU. This imbalance causes\n"
1005 " performance loss. Maybe the domain decomposition limits the PME tuning.\n"
1006 " In that case, try setting the DD grid manually (-dd) or lowering -dds.");
1010 /* We should not end up here, unless the box is
1011 * too small for increasing the cut-off for PME tuning.
1013 GMX_LOG(mdlog
.warning
).asParagraph().appendText(
1014 "NOTE: The GPU has >25% less load than the CPU. This imbalance causes\n"
1015 " performance loss.");
1018 if (gpu_cpu_ratio
> 1.2)
1020 GMX_LOG(mdlog
.warning
).asParagraph().appendText(
1021 "NOTE: The GPU has >20% more load than the CPU. This imbalance causes\n"
1022 " performance loss, consider using a shorter cut-off and a finer PME grid.");
1030 GMX_LOG(mdlog
.warning
).asParagraph().appendText(
1031 "MPI_Barrier was called before each cycle start/stop\n"
1032 "call, so timings are not those of real runs.");
1035 if (wc
->wcc
[ewcNB_XF_BUF_OPS
].n
> 0 &&
1036 (cyc_sum
[ewcDOMDEC
] > tot
*0.1 ||
1037 cyc_sum
[ewcNS
] > tot
*0.1))
1039 /* Only the sim master calls this function, so always print to stderr */
1040 if (wc
->wcc
[ewcDOMDEC
].n
== 0)
1042 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
1043 "NOTE: %d %% of the run time was spent in pair search,\n"
1044 " you might want to increase nstlist (this has no effect on accuracy)\n",
1045 (int)(100*cyc_sum
[ewcNS
]/tot
+0.5));
1049 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
1050 "NOTE: %d %% of the run time was spent in domain decomposition,\n"
1051 " %d %% of the run time was spent in pair search,\n"
1052 " you might want to increase nstlist (this has no effect on accuracy)\n",
1053 (int)(100*cyc_sum
[ewcDOMDEC
]/tot
+0.5),
1054 (int)(100*cyc_sum
[ewcNS
]/tot
+0.5));
1058 if (cyc_sum
[ewcMoveE
] > tot
*0.05)
1060 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
1061 "NOTE: %d %% of the run time was spent communicating energies,\n"
1062 " you might want to use the -gcom option of mdrun\n",
1063 (int)(100*cyc_sum
[ewcMoveE
]/tot
+0.5));
1067 extern gmx_int64_t
wcycle_get_reset_counters(gmx_wallcycle_t wc
)
1074 return wc
->reset_counters
;
1077 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc
, gmx_int64_t reset_counters
)
1084 wc
->reset_counters
= reset_counters
;
1087 void wallcycle_sub_start(gmx_wallcycle_t wc
, int ewcs
)
1089 if (useCycleSubcounters
&& wc
!= nullptr)
1091 wc
->wcsc
[ewcs
].start
= gmx_cycles_read();
1095 void wallcycle_sub_start_nocount(gmx_wallcycle_t wc
, int ewcs
)
1097 if (useCycleSubcounters
&& wc
!= nullptr)
1099 wallcycle_sub_start(wc
, ewcs
);
1104 void wallcycle_sub_stop(gmx_wallcycle_t wc
, int ewcs
)
1106 if (useCycleSubcounters
&& wc
!= nullptr)
1108 wc
->wcsc
[ewcs
].c
+= gmx_cycles_read() - wc
->wcsc
[ewcs
].start
;