3 * This source code is part of
7 * 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-2004, 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 * GROningen Mixture of Alchemy and Childrens' Stories
40 #include "types/commrec.h"
42 #include "gmx_fatal.h"
60 static const t_nrnb_data nbdata[eNRNB] = {
61 /* These are re-used for different NB kernels, since there are so many.
62 * The actual number of flops is set dynamically.
64 { "NB VdW [V&F]", 1 },
66 { "NB Elec. [V&F]", 1 },
67 { "NB Elec. [F]", 1 },
68 { "NB Elec. [W3,V&F]", 1 },
69 { "NB Elec. [W3,F]", 1 },
70 { "NB Elec. [W3-W3,V&F]", 1 },
71 { "NB Elec. [W3-W3,F]", 1 },
72 { "NB Elec. [W4,V&F]", 1 },
73 { "NB Elec. [W4,F]", 1 },
74 { "NB Elec. [W4-W4,V&F]", 1 },
75 { "NB Elec. [W4-W4,F]", 1 },
76 { "NB VdW & Elec. [V&F]", 1 },
77 { "NB VdW & Elec. [F]", 1 },
78 { "NB VdW & Elec. [W3,V&F]", 1 },
79 { "NB VdW & Elec. [W3,F]", 1 },
80 { "NB VdW & Elec. [W3-W3,V&F]", 1 },
81 { "NB VdW & Elec. [W3-W3,F]", 1 },
82 { "NB VdW & Elec. [W4,V&F]", 1 },
83 { "NB VdW & Elec. [W4,F]", 1 },
84 { "NB VdW & Elec. [W4-W4,V&F]", 1 },
85 { "NB VdW & Elec. [W4-W4,F]", 1 },
87 { "NB Generic kernel", 1 },
88 { "NB Generic charge grp kernel", 1 },
89 { "NB Generic AdResS kernel", 1 },
90 { "NB Free energy kernel", 1 },
91 { "NB All-vs-all", 1 },
92 { "NB All-vs-all, GB", 1 },
94 { "Pair Search distance check", 9 }, /* nbnxn pair dist. check */
95 /* nbnxn kernel flops are based on inner-loops without exclusion checks.
96 * Plain Coulomb runs through the RF kernels, except with CUDA.
97 * invsqrt is counted as 6 flops: 1 for _mm_rsqt_ps + 5 for iteration.
98 * The flops are equal for plain-C, x86 SIMD and CUDA, except for:
99 * - plain-C kernel uses one flop more for Coulomb-only (F) than listed
100 * - x86 SIMD LJ geom-comb.rule kernels (fastest) use 2 more flops
101 * - x86 SIMD LJ LB-comb.rule kernels (fast) use 3 (8 for F+E) more flops
102 * - GPU always does exclusions, which requires 2-4 flops, but as invsqrt
103 * is always counted as 6 flops, this roughly compensates.
105 { "NxN RF Elec. + VdW [F]", 38 }, /* nbnxn kernel LJ+RF, no ener */
106 { "NxN RF Elec. + VdW [V&F]", 54 },
107 { "NxN QSTab Elec. + VdW [F]", 41 }, /* nbnxn kernel LJ+tab, no en */
108 { "NxN QSTab Elec. + VdW [V&F]", 59 },
109 { "NxN Ewald Elec. + VdW [F]", 66 }, /* nbnxn kernel LJ+Ewald, no en */
110 { "NxN Ewald Elec. + VdW [V&F]", 107 },
111 { "NxN VdW [F]", 33 }, /* nbnxn kernel LJ, no ener */
112 { "NxN VdW [V&F]", 43 },
113 { "NxN RF Electrostatics [F]", 31 }, /* nbnxn kernel RF, no ener */
114 { "NxN RF Electrostatics [V&F]", 36 },
115 { "NxN QSTab Elec. [F]", 34 }, /* nbnxn kernel tab, no ener */
116 { "NxN QSTab Elec. [V&F]", 41 },
117 { "NxN Ewald Elec. [F]", 61 }, /* nbnxn kernel Ewald, no ener */
118 { "NxN Ewald Elec. [V&F]", 84 },
119 { "1,4 nonbonded interactions", 90 },
120 { "Born radii (Still)", 47 },
121 { "Born radii (HCT/OBC)", 183 },
122 { "Born force chain rule", 15 },
123 { "All-vs-All Still radii", 1 },
124 { "All-vs-All HCT/OBC radii", 1 },
125 { "All-vs-All Born chain rule", 1 },
126 { "Calc Weights", 36 },
128 { "Spread Q Bspline", 2 },
130 { "Gather F Bspline", 6 },
132 { "Convolution", 4 },
135 { "Reset In Box", 3 },
141 { "FENE Bonds", 58 },
142 { "Tab. Bonds", 62 },
143 { "Restraint Potential", 86 },
144 { "Linear Angles", 57 },
146 { "G96Angles", 150 },
147 { "Quartic Angles", 160 },
148 { "Tab. Angles", 169 },
150 { "Impropers", 208 },
151 { "RB-Dihedrals", 247 },
152 { "Four. Dihedrals", 247 },
153 { "Tab. Dihedrals", 227 },
154 { "Dist. Restr.", 200 },
155 { "Orient. Restr.", 200 },
156 { "Dihedral Restr.", 200 },
157 { "Pos. Restr.", 50 },
158 { "Flat-bottom posres", 50 },
159 { "Angle Restr.", 191 },
160 { "Angle Restr. Z", 164 },
161 { "Morse Potent.", 83 },
162 { "Cubic Bonds", 54 },
164 { "Polarization", 59 },
165 { "Anharmonic Polarization", 72 },
166 { "Water Pol.", 62 },
167 { "Thole Pol.", 296 },
170 { "Ext.ens. Update", 54 },
177 { "Constraint-V", 8 },
178 { "Shake-Init", 10 },
179 { "Constraint-Vir", 24 },
181 { "Virtual Site 2", 23 },
182 { "Virtual Site 3", 37 },
183 { "Virtual Site 3fd", 95 },
184 { "Virtual Site 3fad", 176 },
185 { "Virtual Site 3out", 87 },
186 { "Virtual Site 4fd", 110 },
187 { "Virtual Site 4fdn", 254 },
188 { "Virtual Site N", 15 },
189 { "Mixed Generalized Born stuff", 10 }
193 void init_nrnb(t_nrnb *nrnb)
197 for (i = 0; (i < eNRNB); i++)
203 void cp_nrnb(t_nrnb *dest, t_nrnb *src)
207 for (i = 0; (i < eNRNB); i++)
209 dest->n[i] = src->n[i];
213 void add_nrnb(t_nrnb *dest, t_nrnb *s1, t_nrnb *s2)
217 for (i = 0; (i < eNRNB); i++)
219 dest->n[i] = s1->n[i]+s2->n[i];
223 void print_nrnb(FILE *out, t_nrnb *nrnb)
227 for (i = 0; (i < eNRNB); i++)
231 fprintf(out, " %-26s %10.0f.\n", nbdata[i].name, nrnb->n[i]);
236 void _inc_nrnb(t_nrnb *nrnb, int enr, int inc, char *file, int line)
240 printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
241 nbdata[enr].name, enr, inc, file, line);
245 void print_flop(FILE *out, t_nrnb *nrnb, double *nbfs, double *mflop)
248 double mni, frac, tfrac, tflop;
249 const char *myline = "-----------------------------------------------------------------------------";
252 for (i = 0; (i < eNR_NBKERNEL_ALLVSALLGB); i++)
254 if (strstr(nbdata[i].name, "W3-W3") != NULL)
256 *nbfs += 9e-6*nrnb->n[i];
258 else if (strstr(nbdata[i].name, "W3") != NULL)
260 *nbfs += 3e-6*nrnb->n[i];
262 else if (strstr(nbdata[i].name, "W4-W4") != NULL)
264 *nbfs += 10e-6*nrnb->n[i];
266 else if (strstr(nbdata[i].name, "W4") != NULL)
268 *nbfs += 4e-6*nrnb->n[i];
272 *nbfs += 1e-6*nrnb->n[i];
276 for (i = 0; (i < eNRNB); i++)
278 tflop += 1e-6*nrnb->n[i]*nbdata[i].flop;
283 fprintf(out, "No MEGA Flopsen this time\n");
288 fprintf(out, "\n\tM E G A - F L O P S A C C O U N T I N G\n\n");
293 fprintf(out, " NB=Group-cutoff nonbonded kernels NxN=N-by-N cluster Verlet kernels\n");
294 fprintf(out, " RF=Reaction-Field VdW=Van der Waals QSTab=quadratic-spline table\n");
295 fprintf(out, " W3=SPC/TIP3p W4=TIP4p (single or pairs)\n");
296 fprintf(out, " V&F=Potential and force V=Potential only F=Force only\n\n");
298 fprintf(out, " %-32s %16s %15s %7s\n",
299 "Computing:", "M-Number", "M-Flops", "% Flops");
300 fprintf(out, "%s\n", myline);
304 for (i = 0; (i < eNRNB); i++)
306 mni = 1e-6*nrnb->n[i];
307 *mflop += mni*nbdata[i].flop;
308 frac = 100.0*mni*nbdata[i].flop/tflop;
312 fprintf(out, " %-32s %16.6f %15.3f %6.1f\n",
313 nbdata[i].name, mni, mni*nbdata[i].flop, frac);
318 fprintf(out, "%s\n", myline);
319 fprintf(out, " %-32s %16s %15.3f %6.1f\n",
320 "Total", "", *mflop, tfrac);
321 fprintf(out, "%s\n\n", myline);
325 void print_perf(FILE *out, double nodetime, double realtime, int nprocs,
326 gmx_large_int_t nsteps, real delta_t,
327 double nbfs, double mflop,
336 fprintf(out, "%12s %12s %12s %10s\n", "", "Core t (s)", "Wall t (s)", "(%)");
337 fprintf(out, "%12s %12.3f %12.3f %10.1f\n", "Time:",
338 nodetime, realtime, 100.0*nodetime/realtime);
339 /* only print day-hour-sec format if realtime is more than 30 min */
340 if (realtime > 30*60)
342 fprintf(out, "%12s %12s", "", "");
343 pr_difftime(out, realtime);
347 mflop = mflop/realtime;
348 runtime = nsteps*delta_t;
350 if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
352 fprintf(out, "%12s %12s %12s\n",
353 "", "(ns/day)", "(hour/ns)");
354 fprintf(out, "%12s %12.3f %12.3f\n", "Performance:",
355 runtime*24*3.6/realtime, 1000*realtime/(3600*runtime));
359 fprintf(out, "%12s %12s %12s %12s %12s\n",
360 "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
361 "(ns/day)", "(hour/ns)");
362 fprintf(out, "%12s %12.3f %12.3f %12.3f %12.3f\n", "Performance:",
363 nbfs/realtime, (mflop > 1000) ? (mflop/1000) : mflop,
364 runtime*24*3.6/realtime, 1000*realtime/(3600*runtime));
369 if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
371 fprintf(out, "%12s %14s\n",
373 fprintf(out, "%12s %14.1f\n", "Performance:",
374 nsteps*3600.0/realtime);
378 fprintf(out, "%12s %12s %12s %14s\n",
379 "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
381 fprintf(out, "%12s %12.3f %12.3f %14.1f\n", "Performance:",
382 nbfs/realtime, (mflop > 1000) ? (mflop/1000) : mflop,
383 nsteps*3600.0/realtime);
389 int cost_nrnb(int enr)
391 return nbdata[enr].flop;
394 const char *nrnb_str(int enr)
396 return nbdata[enr].name;
399 static const int force_index[] = {
400 eNR_BONDS, eNR_ANGLES, eNR_PROPER, eNR_IMPROPER,
401 eNR_RB, eNR_DISRES, eNR_ORIRES, eNR_POSRES,
402 eNR_FBPOSRES, eNR_NS,
404 #define NFORCE_INDEX asize(force_index)
406 static const int constr_index[] = {
407 eNR_SHAKE, eNR_SHAKE_RIJ, eNR_SETTLE, eNR_UPDATE, eNR_PCOUPL,
408 eNR_CONSTR_VIR, eNR_CONSTR_V
410 #define NCONSTR_INDEX asize(constr_index)
412 static double pr_av(FILE *log, t_commrec *cr,
413 double fav, double ftot[], const char *title)
421 fav /= cr->nnodes - cr->npmenodes;
422 fprintf(log, "\n %-26s", title);
423 for (i = 0; (i < cr->nnodes); i++)
425 dperc = (100.0*ftot[i])/fav;
426 unb = max(unb, dperc);
428 fprintf(log, "%3d ", perc);
433 fprintf(log, "%6d%%\n\n", perc);
437 fprintf(log, "\n\n");
443 void pr_load(FILE *log, t_commrec *cr, t_nrnb nrnb[])
446 double dperc, unb, uf, us;
452 snew(ftot, cr->nnodes);
453 snew(stot, cr->nnodes);
455 for (i = 0; (i < cr->nnodes); i++)
457 add_nrnb(av, av, &(nrnb[i]));
458 /* Cost due to forces */
459 for (j = 0; (j < eNR_NBKERNEL_ALLVSALLGB); j++)
461 ftot[i] += nrnb[i].n[j]*cost_nrnb(j);
463 for (j = 0; (j < NFORCE_INDEX); j++)
465 ftot[i] += nrnb[i].n[force_index[j]]*cost_nrnb(force_index[j]);
468 for (j = 0; (j < NCONSTR_INDEX); j++)
470 stot[i] += nrnb[i].n[constr_index[j]]*cost_nrnb(constr_index[j]);
473 for (j = 0; (j < eNRNB); j++)
475 av->n[j] = av->n[j]/(double)(cr->nnodes - cr->npmenodes);
478 fprintf(log, "\nDetailed load balancing info in percentage of average\n");
480 fprintf(log, " Type NODE:");
481 for (i = 0; (i < cr->nnodes); i++)
483 fprintf(log, "%3d ", i);
485 fprintf(log, "Scaling\n");
486 fprintf(log, "---------------------------");
487 for (i = 0; (i < cr->nnodes); i++)
489 fprintf(log, "----");
491 fprintf(log, "-------\n");
493 for (j = 0; (j < eNRNB); j++)
498 fprintf(log, " %-26s", nrnb_str(j));
499 for (i = 0; (i < cr->nnodes); i++)
501 dperc = (100.0*nrnb[i].n[j])/av->n[j];
502 unb = max(unb, dperc);
504 fprintf(log, "%3d ", perc);
509 fprintf(log, "%6d%%\n", perc);
518 for (i = 0; (i < cr->nnodes); i++)
523 uf = pr_av(log, cr, fav, ftot, "Total Force");
524 us = pr_av(log, cr, sav, stot, "Total Constr.");
526 unb = (uf*fav+us*sav)/(fav+sav);
530 fprintf(log, "\nTotal Scaling: %.0f%% of max performance\n\n", unb);