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-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2017,2018 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source 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.
48 #include "gromacs/mdtypes/commrec.h"
49 #include "gromacs/mdtypes/md_enums.h"
50 #include "gromacs/utility/arraysize.h"
59 static const t_nrnb_data nbdata[eNRNB] = {
60 /* These are re-used for different NB kernels, since there are so many.
61 * The actual number of flops is set dynamically.
63 { "NB VdW [V&F]", 1 },
65 { "NB Elec. [V&F]", 1 },
66 { "NB Elec. [F]", 1 },
67 { "NB Elec. [W3,V&F]", 1 },
68 { "NB Elec. [W3,F]", 1 },
69 { "NB Elec. [W3-W3,V&F]", 1 },
70 { "NB Elec. [W3-W3,F]", 1 },
71 { "NB Elec. [W4,V&F]", 1 },
72 { "NB Elec. [W4,F]", 1 },
73 { "NB Elec. [W4-W4,V&F]", 1 },
74 { "NB Elec. [W4-W4,F]", 1 },
75 { "NB VdW & Elec. [V&F]", 1 },
76 { "NB VdW & Elec. [F]", 1 },
77 { "NB VdW & Elec. [W3,V&F]", 1 },
78 { "NB VdW & Elec. [W3,F]", 1 },
79 { "NB VdW & Elec. [W3-W3,V&F]", 1 },
80 { "NB VdW & Elec. [W3-W3,F]", 1 },
81 { "NB VdW & Elec. [W4,V&F]", 1 },
82 { "NB VdW & Elec. [W4,F]", 1 },
83 { "NB VdW & Elec. [W4-W4,V&F]", 1 },
84 { "NB VdW & Elec. [W4-W4,F]", 1 },
86 { "NB Generic kernel", 1 },
87 { "NB Generic charge grp kernel", 1 },
88 { "NB Free energy kernel", 1 },
90 { "Pair Search distance check", 9 }, /* nbnxn pair dist. check */
91 /* nbnxn kernel flops are based on inner-loops without exclusion checks.
92 * Plain Coulomb runs through the RF kernels, except with GPUs.
93 * invsqrt is counted as 6 flops: 1 for _mm_rsqt_ps + 5 for iteration.
94 * The flops are equal for plain-C, x86 SIMD and GPUs, except for:
95 * - plain-C kernel uses one flop more for Coulomb-only (F) than listed
96 * - x86 SIMD LJ geom-comb.rule kernels (fastest) use 2 more flops
97 * - x86 SIMD LJ LB-comb.rule kernels (fast) use 3 (8 for F+E) more flops
98 * - GPU always does exclusions, which requires 2-4 flops, but as invsqrt
99 * is always counted as 6 flops, this roughly compensates.
101 { "NxN RF Elec. + LJ [F]", 38 }, /* nbnxn kernel LJ+RF, no ener */
102 { "NxN RF Elec. + LJ [V&F]", 54 },
103 { "NxN QSTab Elec. + LJ [F]", 41 }, /* nbnxn kernel LJ+tab, no en */
104 { "NxN QSTab Elec. + LJ [V&F]", 59 },
105 { "NxN Ewald Elec. + LJ [F]", 66 }, /* nbnxn kernel LJ+Ewald, no en */
106 { "NxN Ewald Elec. + LJ [V&F]", 107 },
107 { "NxN LJ [F]", 33 }, /* nbnxn kernel LJ, no ener */
108 { "NxN LJ [V&F]", 43 },
109 { "NxN RF Electrostatics [F]", 31 }, /* nbnxn kernel RF, no ener */
110 { "NxN RF Electrostatics [V&F]", 36 },
111 { "NxN QSTab Elec. [F]", 34 }, /* nbnxn kernel tab, no ener */
112 { "NxN QSTab Elec. [V&F]", 41 },
113 { "NxN Ewald Elec. [F]", 61 }, /* nbnxn kernel Ewald, no ener */
114 { "NxN Ewald Elec. [V&F]", 84 },
115 /* The switch function flops should be added to the LJ kernels above */
116 { "NxN LJ add F-switch [F]", 12 }, /* extra cost for LJ F-switch */
117 { "NxN LJ add F-switch [V&F]", 22 },
118 { "NxN LJ add P-switch [F]", 27 }, /* extra cost for LJ P-switch */
119 { "NxN LJ add P-switch [V&F]", 20 },
120 { "NxN LJ add LJ Ewald [F]", 36 }, /* extra cost for LJ Ewald */
121 { "NxN LJ add LJ Ewald [V&F]", 33 },
122 { "1,4 nonbonded interactions", 90 },
123 { "Calc Weights", 36 },
125 { "Spread Q Bspline", 2 },
127 { "Gather F Bspline", 6 },
129 { "Convolution", 4 },
132 { "Reset In Box", 3 },
138 { "FENE Bonds", 58 },
139 { "Tab. Bonds", 62 },
140 { "Restraint Potential", 86 },
141 { "Linear Angles", 57 },
143 { "G96Angles", 150 },
144 { "Quartic Angles", 160 },
145 { "Tab. Angles", 169 },
147 { "Impropers", 208 },
148 { "RB-Dihedrals", 247 },
149 { "Four. Dihedrals", 247 },
150 { "Tab. Dihedrals", 227 },
151 { "Dist. Restr.", 200 },
152 { "Orient. Restr.", 200 },
153 { "Dihedral Restr.", 200 },
154 { "Pos. Restr.", 50 },
155 { "Flat-bottom posres", 50 },
156 { "Angle Restr.", 191 },
157 { "Angle Restr. Z", 164 },
158 { "Morse Potent.", 83 },
159 { "Cubic Bonds", 54 },
161 { "Polarization", 59 },
162 { "Anharmonic Polarization", 72 },
163 { "Water Pol.", 62 },
164 { "Thole Pol.", 296 },
167 { "Ext.ens. Update", 54 },
174 { "Constraint-V", 9 },
175 { "Shake-Init", 10 },
176 { "Constraint-Vir", 24 },
178 { "Virtual Site 2", 23 },
179 { "Virtual Site 2fd", 63 },
180 { "Virtual Site 3", 37 },
181 { "Virtual Site 3fd", 95 },
182 { "Virtual Site 3fad", 176 },
183 { "Virtual Site 3out", 87 },
184 { "Virtual Site 4fd", 110 },
185 { "Virtual Site 4fdn", 254 },
186 { "Virtual Site N", 15 },
187 { "CMAP", 1700 }, // Estimate!
188 { "Urey-Bradley", 183 },
189 { "Cross-Bond-Bond", 163 },
190 { "Cross-Bond-Angle", 163 }
193 static void pr_two(FILE* out, int c, int i)
197 fprintf(out, "%c0%1d", c, i);
201 fprintf(out, "%c%2d", c, i);
205 static void pr_difftime(FILE* out, double dt)
207 int ndays, nhours, nmins, nsecs;
208 gmx_bool bPrint, bPrinted;
210 ndays = static_cast<int>(dt / (24 * 3600));
211 dt = dt - 24 * 3600 * ndays;
212 nhours = static_cast<int>(dt / 3600);
213 dt = dt - 3600 * nhours;
214 nmins = static_cast<int>(dt / 60);
215 dt = dt - nmins * 60;
216 nsecs = static_cast<int>(dt);
217 bPrint = (ndays > 0);
221 fprintf(out, "%d", ndays);
223 bPrint = bPrint || (nhours > 0);
228 pr_two(out, 'd', nhours);
232 fprintf(out, "%d", nhours);
235 bPrinted = bPrinted || bPrint;
236 bPrint = bPrint || (nmins > 0);
241 pr_two(out, 'h', nmins);
245 fprintf(out, "%d", nmins);
248 bPrinted = bPrinted || bPrint;
251 pr_two(out, ':', nsecs);
255 fprintf(out, "%ds", nsecs);
260 void clear_nrnb(t_nrnb* nrnb)
264 for (i = 0; (i < eNRNB); i++)
270 void add_nrnb(t_nrnb* dest, t_nrnb* s1, t_nrnb* s2)
274 for (i = 0; (i < eNRNB); i++)
276 dest->n[i] = s1->n[i] + s2->n[i];
280 void print_nrnb(FILE* out, t_nrnb* nrnb)
284 for (i = 0; (i < eNRNB); i++)
288 fprintf(out, " %-26s %10.0f.\n", nbdata[i].name, nrnb->n[i]);
293 void _inc_nrnb(t_nrnb* nrnb, int enr, int inc, char gmx_unused* file, int gmx_unused line)
297 printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n", nbdata[enr].name, enr, inc,
302 /* Returns in enr is the index of a full nbnxn VdW kernel */
303 static gmx_bool nrnb_is_nbnxn_vdw_kernel(int enr)
305 return (enr >= eNR_NBNXN_LJ_RF && enr <= eNR_NBNXN_LJ_E);
308 /* Returns in enr is the index of an nbnxn kernel addition (LJ modification) */
309 static gmx_bool nrnb_is_nbnxn_kernel_addition(int enr)
311 return (enr >= eNR_NBNXN_ADD_LJ_FSW && enr <= eNR_NBNXN_ADD_LJ_EWALD_E);
314 void print_flop(FILE* out, t_nrnb* nrnb, double* nbfs, double* mflop)
317 double mni, frac, tfrac, tflop;
319 "-----------------------------------------------------------------------------";
322 for (i = 0; (i < eNR_NBKERNEL_TOTAL_NR); i++)
324 if (std::strstr(nbdata[i].name, "W3-W3") != nullptr)
326 *nbfs += 9e-6 * nrnb->n[i];
328 else if (std::strstr(nbdata[i].name, "W3") != nullptr)
330 *nbfs += 3e-6 * nrnb->n[i];
332 else if (std::strstr(nbdata[i].name, "W4-W4") != nullptr)
334 *nbfs += 10e-6 * nrnb->n[i];
336 else if (std::strstr(nbdata[i].name, "W4") != nullptr)
338 *nbfs += 4e-6 * nrnb->n[i];
342 *nbfs += 1e-6 * nrnb->n[i];
346 for (i = 0; (i < eNRNB); i++)
348 tflop += 1e-6 * nrnb->n[i] * nbdata[i].flop;
353 fprintf(out, "No MEGA Flopsen this time\n");
358 fprintf(out, "\n\tM E G A - F L O P S A C C O U N T I N G\n\n");
363 fprintf(out, " NB=Group-cutoff nonbonded kernels NxN=N-by-N cluster Verlet kernels\n");
364 fprintf(out, " RF=Reaction-Field VdW=Van der Waals QSTab=quadratic-spline table\n");
365 fprintf(out, " W3=SPC/TIP3p W4=TIP4p (single or pairs)\n");
366 fprintf(out, " V&F=Potential and force V=Potential only F=Force only\n\n");
368 fprintf(out, " %-32s %16s %15s %7s\n", "Computing:", "M-Number", "M-Flops", "% Flops");
369 fprintf(out, "%s\n", myline);
373 for (i = 0; (i < eNRNB); i++)
375 mni = 1e-6 * nrnb->n[i];
376 /* Skip empty entries and nbnxn additional flops,
377 * which have been added to the kernel entry.
379 if (mni > 0 && !nrnb_is_nbnxn_kernel_addition(i))
383 flop = nbdata[i].flop;
384 if (nrnb_is_nbnxn_vdw_kernel(i))
386 /* Possibly add the cost of an LJ switch/Ewald function */
387 for (j = eNR_NBNXN_ADD_LJ_FSW; j <= eNR_NBNXN_ADD_LJ_EWALD; j += 2)
391 /* Select the force or energy flop count */
392 e_kernel_add = j + ((i - eNR_NBNXN_LJ_RF) % 2);
394 if (nrnb->n[e_kernel_add] > 0)
396 flop += nbdata[e_kernel_add].flop;
400 *mflop += mni * flop;
401 frac = 100.0 * mni * flop / tflop;
405 fprintf(out, " %-32s %16.6f %15.3f %6.1f\n", nbdata[i].name, mni, mni * flop, frac);
411 fprintf(out, "%s\n", myline);
412 fprintf(out, " %-32s %16s %15.3f %6.1f\n", "Total", "", *mflop, tfrac);
413 fprintf(out, "%s\n\n", myline);
415 if (nrnb->n[eNR_NBKERNEL_GENERIC] > 0)
418 "WARNING: Using the slow generic C kernel. This is fine if you are\n"
419 "comparing different implementations or MD software. Routine\n"
420 "simulations should use a different non-bonded setup for much better\n"
426 void print_perf(FILE* out,
427 double time_per_thread,
428 double time_per_node,
434 double wallclocktime;
438 if (time_per_node > 0)
440 fprintf(out, "%12s %12s %12s %10s\n", "", "Core t (s)", "Wall t (s)", "(%)");
441 fprintf(out, "%12s %12.3f %12.3f %10.1f\n", "Time:", time_per_thread, time_per_node,
442 100.0 * time_per_thread / time_per_node);
443 /* only print day-hour-sec format if time_per_node is more than 30 min */
444 if (time_per_node > 30 * 60)
446 fprintf(out, "%12s %12s", "", "");
447 pr_difftime(out, time_per_node);
451 mflop = mflop / time_per_node;
452 wallclocktime = nsteps * delta_t;
454 if (getenv("GMX_DETAILED_PERF_STATS") == nullptr)
456 fprintf(out, "%12s %12s %12s\n", "", "(ns/day)", "(hour/ns)");
457 fprintf(out, "%12s %12.3f %12.3f\n", "Performance:", wallclocktime * 24 * 3.6 / time_per_node,
458 1000 * time_per_node / (3600 * wallclocktime));
462 fprintf(out, "%12s %12s %12s %12s %12s\n", "", "(Mnbf/s)",
463 (mflop > 1000) ? "(GFlops)" : "(MFlops)", "(ns/day)", "(hour/ns)");
464 fprintf(out, "%12s %12.3f %12.3f %12.3f %12.3f\n", "Performance:", nbfs / time_per_node,
465 (mflop > 1000) ? (mflop / 1000) : mflop, wallclocktime * 24 * 3.6 / time_per_node,
466 1000 * time_per_node / (3600 * wallclocktime));
471 if (getenv("GMX_DETAILED_PERF_STATS") == nullptr)
473 fprintf(out, "%12s %14s\n", "", "(steps/hour)");
474 fprintf(out, "%12s %14.1f\n", "Performance:", nsteps * 3600.0 / time_per_node);
478 fprintf(out, "%12s %12s %12s %14s\n", "", "(Mnbf/s)",
479 (mflop > 1000) ? "(GFlops)" : "(MFlops)", "(steps/hour)");
480 fprintf(out, "%12s %12.3f %12.3f %14.1f\n", "Performance:", nbfs / time_per_node,
481 (mflop > 1000) ? (mflop / 1000) : mflop, nsteps * 3600.0 / time_per_node);
487 int cost_nrnb(int enr)
489 return nbdata[enr].flop;
492 const char* nrnb_str(int enr)
494 return nbdata[enr].name;
497 static const int force_index[] = {
498 eNR_BONDS, eNR_ANGLES, eNR_PROPER, eNR_IMPROPER, eNR_RB,
499 eNR_DISRES, eNR_ORIRES, eNR_POSRES, eNR_FBPOSRES, eNR_NS,
501 #define NFORCE_INDEX asize(force_index)
503 static const int constr_index[] = { eNR_SHAKE, eNR_SHAKE_RIJ, eNR_SETTLE, eNR_UPDATE,
504 eNR_PCOUPL, eNR_CONSTR_VIR, eNR_CONSTR_V };
505 #define NCONSTR_INDEX asize(constr_index)
507 static double pr_av(FILE* log, t_commrec* cr, double fav, const double ftot[], const char* title)
515 fav /= cr->nnodes - cr->npmenodes;
516 fprintf(log, "\n %-26s", title);
517 for (i = 0; (i < cr->nnodes); i++)
519 dperc = (100.0 * ftot[i]) / fav;
520 unb = std::max(unb, dperc);
521 perc = static_cast<int>(dperc);
522 fprintf(log, "%3d ", perc);
526 perc = static_cast<int>(10000.0 / unb);
527 fprintf(log, "%6d%%\n\n", perc);
531 fprintf(log, "\n\n");
537 void pr_load(FILE* log, t_commrec* cr, t_nrnb nrnb[])
541 std::vector<double> ftot(cr->nnodes);
542 std::vector<double> stot(cr->nnodes);
543 for (int i = 0; (i < cr->nnodes); i++)
545 add_nrnb(&av, &av, &(nrnb[i]));
546 /* Cost due to forces */
547 for (int j = 0; (j < eNR_NBKERNEL_TOTAL_NR); j++)
549 ftot[i] += nrnb[i].n[j] * cost_nrnb(j);
551 for (int j = 0; (j < NFORCE_INDEX); j++)
553 ftot[i] += nrnb[i].n[force_index[j]] * cost_nrnb(force_index[j]);
556 for (int j = 0; (j < NCONSTR_INDEX); j++)
558 stot[i] += nrnb[i].n[constr_index[j]] * cost_nrnb(constr_index[j]);
561 for (int j = 0; (j < eNRNB); j++)
563 av.n[j] = av.n[j] / static_cast<double>(cr->nnodes - cr->npmenodes);
566 fprintf(log, "\nDetailed load balancing info in percentage of average\n");
568 fprintf(log, " Type RANK:");
569 for (int i = 0; (i < cr->nnodes); i++)
571 fprintf(log, "%3d ", i);
573 fprintf(log, "Scaling\n");
574 fprintf(log, "---------------------------");
575 for (int i = 0; (i < cr->nnodes); i++)
577 fprintf(log, "----");
579 fprintf(log, "-------\n");
581 for (int j = 0; (j < eNRNB); j++)
588 fprintf(log, " %-26s", nrnb_str(j));
589 for (int i = 0; (i < cr->nnodes); i++)
591 dperc = (100.0 * nrnb[i].n[j]) / av.n[j];
592 unb = std::max(unb, dperc);
593 perc = static_cast<int>(dperc);
594 fprintf(log, "%3d ", perc);
598 perc = static_cast<int>(10000.0 / unb);
599 fprintf(log, "%6d%%\n", perc);
609 for (int i = 0; (i < cr->nnodes); i++)
614 double uf = pr_av(log, cr, fav, ftot.data(), "Total Force");
615 double us = pr_av(log, cr, sav, stot.data(), "Total Constr.");
617 double unb = (uf * fav + us * sav) / (fav + sav);
621 fprintf(log, "\nTotal Scaling: %.0f%% of max performance\n\n", unb);