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, 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.
42 #include "types/commrec.h"
44 #include "gmx_fatal.h"
57 static const t_nrnb_data nbdata[eNRNB] = {
58 /* These are re-used for different NB kernels, since there are so many.
59 * The actual number of flops is set dynamically.
61 { "NB VdW [V&F]", 1 },
63 { "NB Elec. [V&F]", 1 },
64 { "NB Elec. [F]", 1 },
65 { "NB Elec. [W3,V&F]", 1 },
66 { "NB Elec. [W3,F]", 1 },
67 { "NB Elec. [W3-W3,V&F]", 1 },
68 { "NB Elec. [W3-W3,F]", 1 },
69 { "NB Elec. [W4,V&F]", 1 },
70 { "NB Elec. [W4,F]", 1 },
71 { "NB Elec. [W4-W4,V&F]", 1 },
72 { "NB Elec. [W4-W4,F]", 1 },
73 { "NB VdW & Elec. [V&F]", 1 },
74 { "NB VdW & Elec. [F]", 1 },
75 { "NB VdW & Elec. [W3,V&F]", 1 },
76 { "NB VdW & Elec. [W3,F]", 1 },
77 { "NB VdW & Elec. [W3-W3,V&F]", 1 },
78 { "NB VdW & Elec. [W3-W3,F]", 1 },
79 { "NB VdW & Elec. [W4,V&F]", 1 },
80 { "NB VdW & Elec. [W4,F]", 1 },
81 { "NB VdW & Elec. [W4-W4,V&F]", 1 },
82 { "NB VdW & Elec. [W4-W4,F]", 1 },
84 { "NB Generic kernel", 1 },
85 { "NB Generic charge grp kernel", 1 },
86 { "NB Generic AdResS kernel", 1 },
87 { "NB Free energy kernel", 1 },
88 { "NB All-vs-all", 1 },
89 { "NB All-vs-all, GB", 1 },
91 { "Pair Search distance check", 9 }, /* nbnxn pair dist. check */
92 /* nbnxn kernel flops are based on inner-loops without exclusion checks.
93 * Plain Coulomb runs through the RF kernels, except with CUDA.
94 * invsqrt is counted as 6 flops: 1 for _mm_rsqt_ps + 5 for iteration.
95 * The flops are equal for plain-C, x86 SIMD and CUDA, except for:
96 * - plain-C kernel uses one flop more for Coulomb-only (F) than listed
97 * - x86 SIMD LJ geom-comb.rule kernels (fastest) use 2 more flops
98 * - x86 SIMD LJ LB-comb.rule kernels (fast) use 3 (8 for F+E) more flops
99 * - GPU always does exclusions, which requires 2-4 flops, but as invsqrt
100 * is always counted as 6 flops, this roughly compensates.
102 { "NxN RF Elec. + VdW [F]", 38 }, /* nbnxn kernel LJ+RF, no ener */
103 { "NxN RF Elec. + VdW [V&F]", 54 },
104 { "NxN QSTab Elec. + VdW [F]", 41 }, /* nbnxn kernel LJ+tab, no en */
105 { "NxN QSTab Elec. + VdW [V&F]", 59 },
106 { "NxN Ewald Elec. + VdW [F]", 66 }, /* nbnxn kernel LJ+Ewald, no en */
107 { "NxN Ewald Elec. + VdW [V&F]", 107 },
108 { "NxN VdW [F]", 33 }, /* nbnxn kernel LJ, no ener */
109 { "NxN VdW [V&F]", 43 },
110 { "NxN RF Electrostatics [F]", 31 }, /* nbnxn kernel RF, no ener */
111 { "NxN RF Electrostatics [V&F]", 36 },
112 { "NxN QSTab Elec. [F]", 34 }, /* nbnxn kernel tab, no ener */
113 { "NxN QSTab Elec. [V&F]", 41 },
114 { "NxN Ewald Elec. [F]", 61 }, /* nbnxn kernel Ewald, no ener */
115 { "NxN Ewald Elec. [V&F]", 84 },
116 { "1,4 nonbonded interactions", 90 },
117 { "Born radii (Still)", 47 },
118 { "Born radii (HCT/OBC)", 183 },
119 { "Born force chain rule", 15 },
120 { "All-vs-All Still radii", 1 },
121 { "All-vs-All HCT/OBC radii", 1 },
122 { "All-vs-All Born chain rule", 1 },
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", 8 },
175 { "Shake-Init", 10 },
176 { "Constraint-Vir", 24 },
178 { "Virtual Site 2", 23 },
179 { "Virtual Site 3", 37 },
180 { "Virtual Site 3fd", 95 },
181 { "Virtual Site 3fad", 176 },
182 { "Virtual Site 3out", 87 },
183 { "Virtual Site 4fd", 110 },
184 { "Virtual Site 4fdn", 254 },
185 { "Virtual Site N", 15 },
186 { "Mixed Generalized Born stuff", 10 }
189 static void pr_two(FILE *out, int c, int i)
193 fprintf(out, "%c0%1d", c, i);
197 fprintf(out, "%c%2d", c, i);
201 static void pr_difftime(FILE *out, double dt)
203 int ndays, nhours, nmins, nsecs;
204 gmx_bool bPrint, bPrinted;
206 ndays = dt/(24*3600);
207 dt = dt-24*3600*ndays;
213 bPrint = (ndays > 0);
217 fprintf(out, "%d", ndays);
219 bPrint = bPrint || (nhours > 0);
224 pr_two(out, 'd', nhours);
228 fprintf(out, "%d", nhours);
231 bPrinted = bPrinted || bPrint;
232 bPrint = bPrint || (nmins > 0);
237 pr_two(out, 'h', nmins);
241 fprintf(out, "%d", nmins);
244 bPrinted = bPrinted || bPrint;
247 pr_two(out, ':', nsecs);
251 fprintf(out, "%ds", nsecs);
256 void init_nrnb(t_nrnb *nrnb)
260 for (i = 0; (i < eNRNB); i++)
266 void cp_nrnb(t_nrnb *dest, t_nrnb *src)
270 for (i = 0; (i < eNRNB); i++)
272 dest->n[i] = src->n[i];
276 void add_nrnb(t_nrnb *dest, t_nrnb *s1, t_nrnb *s2)
280 for (i = 0; (i < eNRNB); i++)
282 dest->n[i] = s1->n[i]+s2->n[i];
286 void print_nrnb(FILE *out, t_nrnb *nrnb)
290 for (i = 0; (i < eNRNB); i++)
294 fprintf(out, " %-26s %10.0f.\n", nbdata[i].name, nrnb->n[i]);
299 void _inc_nrnb(t_nrnb *nrnb, int enr, int inc, char gmx_unused *file, int gmx_unused line)
303 printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
304 nbdata[enr].name, enr, inc, file, line);
308 void print_flop(FILE *out, t_nrnb *nrnb, double *nbfs, double *mflop)
311 double mni, frac, tfrac, tflop;
312 const char *myline = "-----------------------------------------------------------------------------";
315 for (i = 0; (i < eNR_NBKERNEL_ALLVSALLGB); i++)
317 if (strstr(nbdata[i].name, "W3-W3") != NULL)
319 *nbfs += 9e-6*nrnb->n[i];
321 else if (strstr(nbdata[i].name, "W3") != NULL)
323 *nbfs += 3e-6*nrnb->n[i];
325 else if (strstr(nbdata[i].name, "W4-W4") != NULL)
327 *nbfs += 10e-6*nrnb->n[i];
329 else if (strstr(nbdata[i].name, "W4") != NULL)
331 *nbfs += 4e-6*nrnb->n[i];
335 *nbfs += 1e-6*nrnb->n[i];
339 for (i = 0; (i < eNRNB); i++)
341 tflop += 1e-6*nrnb->n[i]*nbdata[i].flop;
346 fprintf(out, "No MEGA Flopsen this time\n");
351 fprintf(out, "\n\tM E G A - F L O P S A C C O U N T I N G\n\n");
356 fprintf(out, " NB=Group-cutoff nonbonded kernels NxN=N-by-N cluster Verlet kernels\n");
357 fprintf(out, " RF=Reaction-Field VdW=Van der Waals QSTab=quadratic-spline table\n");
358 fprintf(out, " W3=SPC/TIP3p W4=TIP4p (single or pairs)\n");
359 fprintf(out, " V&F=Potential and force V=Potential only F=Force only\n\n");
361 fprintf(out, " %-32s %16s %15s %7s\n",
362 "Computing:", "M-Number", "M-Flops", "% Flops");
363 fprintf(out, "%s\n", myline);
367 for (i = 0; (i < eNRNB); i++)
369 mni = 1e-6*nrnb->n[i];
370 *mflop += mni*nbdata[i].flop;
371 frac = 100.0*mni*nbdata[i].flop/tflop;
375 fprintf(out, " %-32s %16.6f %15.3f %6.1f\n",
376 nbdata[i].name, mni, mni*nbdata[i].flop, frac);
381 fprintf(out, "%s\n", myline);
382 fprintf(out, " %-32s %16s %15.3f %6.1f\n",
383 "Total", "", *mflop, tfrac);
384 fprintf(out, "%s\n\n", myline);
388 void print_perf(FILE *out, double time_per_thread, double time_per_node,
389 gmx_int64_t nsteps, real delta_t,
390 double nbfs, double mflop)
396 if (time_per_node > 0)
398 fprintf(out, "%12s %12s %12s %10s\n", "", "Core t (s)", "Wall t (s)", "(%)");
399 fprintf(out, "%12s %12.3f %12.3f %10.1f\n", "Time:",
400 time_per_thread, time_per_node, 100.0*time_per_thread/time_per_node);
401 /* only print day-hour-sec format if time_per_node is more than 30 min */
402 if (time_per_node > 30*60)
404 fprintf(out, "%12s %12s", "", "");
405 pr_difftime(out, time_per_node);
409 mflop = mflop/time_per_node;
410 wallclocktime = nsteps*delta_t;
412 if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
414 fprintf(out, "%12s %12s %12s\n",
415 "", "(ns/day)", "(hour/ns)");
416 fprintf(out, "%12s %12.3f %12.3f\n", "Performance:",
417 wallclocktime*24*3.6/time_per_node, 1000*time_per_node/(3600*wallclocktime));
421 fprintf(out, "%12s %12s %12s %12s %12s\n",
422 "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
423 "(ns/day)", "(hour/ns)");
424 fprintf(out, "%12s %12.3f %12.3f %12.3f %12.3f\n", "Performance:",
425 nbfs/time_per_node, (mflop > 1000) ? (mflop/1000) : mflop,
426 wallclocktime*24*3.6/time_per_node, 1000*time_per_node/(3600*wallclocktime));
431 if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
433 fprintf(out, "%12s %14s\n",
435 fprintf(out, "%12s %14.1f\n", "Performance:",
436 nsteps*3600.0/time_per_node);
440 fprintf(out, "%12s %12s %12s %14s\n",
441 "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
443 fprintf(out, "%12s %12.3f %12.3f %14.1f\n", "Performance:",
444 nbfs/time_per_node, (mflop > 1000) ? (mflop/1000) : mflop,
445 nsteps*3600.0/time_per_node);
451 int cost_nrnb(int enr)
453 return nbdata[enr].flop;
456 const char *nrnb_str(int enr)
458 return nbdata[enr].name;
461 static const int force_index[] = {
462 eNR_BONDS, eNR_ANGLES, eNR_PROPER, eNR_IMPROPER,
463 eNR_RB, eNR_DISRES, eNR_ORIRES, eNR_POSRES,
464 eNR_FBPOSRES, eNR_NS,
466 #define NFORCE_INDEX asize(force_index)
468 static const int constr_index[] = {
469 eNR_SHAKE, eNR_SHAKE_RIJ, eNR_SETTLE, eNR_UPDATE, eNR_PCOUPL,
470 eNR_CONSTR_VIR, eNR_CONSTR_V
472 #define NCONSTR_INDEX asize(constr_index)
474 static double pr_av(FILE *log, t_commrec *cr,
475 double fav, double ftot[], const char *title)
483 fav /= cr->nnodes - cr->npmenodes;
484 fprintf(log, "\n %-26s", title);
485 for (i = 0; (i < cr->nnodes); i++)
487 dperc = (100.0*ftot[i])/fav;
488 unb = max(unb, dperc);
490 fprintf(log, "%3d ", perc);
495 fprintf(log, "%6d%%\n\n", perc);
499 fprintf(log, "\n\n");
505 void pr_load(FILE *log, t_commrec *cr, t_nrnb nrnb[])
508 double dperc, unb, uf, us;
514 snew(ftot, cr->nnodes);
515 snew(stot, cr->nnodes);
517 for (i = 0; (i < cr->nnodes); i++)
519 add_nrnb(av, av, &(nrnb[i]));
520 /* Cost due to forces */
521 for (j = 0; (j < eNR_NBKERNEL_ALLVSALLGB); j++)
523 ftot[i] += nrnb[i].n[j]*cost_nrnb(j);
525 for (j = 0; (j < NFORCE_INDEX); j++)
527 ftot[i] += nrnb[i].n[force_index[j]]*cost_nrnb(force_index[j]);
530 for (j = 0; (j < NCONSTR_INDEX); j++)
532 stot[i] += nrnb[i].n[constr_index[j]]*cost_nrnb(constr_index[j]);
535 for (j = 0; (j < eNRNB); j++)
537 av->n[j] = av->n[j]/(double)(cr->nnodes - cr->npmenodes);
540 fprintf(log, "\nDetailed load balancing info in percentage of average\n");
542 fprintf(log, " Type NODE:");
543 for (i = 0; (i < cr->nnodes); i++)
545 fprintf(log, "%3d ", i);
547 fprintf(log, "Scaling\n");
548 fprintf(log, "---------------------------");
549 for (i = 0; (i < cr->nnodes); i++)
551 fprintf(log, "----");
553 fprintf(log, "-------\n");
555 for (j = 0; (j < eNRNB); j++)
560 fprintf(log, " %-26s", nrnb_str(j));
561 for (i = 0; (i < cr->nnodes); i++)
563 dperc = (100.0*nrnb[i].n[j])/av->n[j];
564 unb = max(unb, dperc);
566 fprintf(log, "%3d ", perc);
571 fprintf(log, "%6d%%\n", perc);
580 for (i = 0; (i < cr->nnodes); i++)
585 uf = pr_av(log, cr, fav, ftot, "Total Force");
586 us = pr_av(log, cr, sav, stot, "Total Constr.");
588 unb = (uf*fav+us*sav)/(fav+sav);
592 fprintf(log, "\nTotal Scaling: %.0f%% of max performance\n\n", unb);