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, 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 "gromacs/legacyheaders/nrnb.h"
46 #include "gromacs/legacyheaders/names.h"
47 #include "gromacs/legacyheaders/types/commrec.h"
48 #include "gromacs/legacyheaders/types/nrnb.h"
49 #include "gromacs/utility/arraysize.h"
50 #include "gromacs/utility/smalloc.h"
58 static const t_nrnb_data nbdata
[eNRNB
] = {
59 /* These are re-used for different NB kernels, since there are so many.
60 * The actual number of flops is set dynamically.
62 { "NB VdW [V&F]", 1 },
64 { "NB Elec. [V&F]", 1 },
65 { "NB Elec. [F]", 1 },
66 { "NB Elec. [W3,V&F]", 1 },
67 { "NB Elec. [W3,F]", 1 },
68 { "NB Elec. [W3-W3,V&F]", 1 },
69 { "NB Elec. [W3-W3,F]", 1 },
70 { "NB Elec. [W4,V&F]", 1 },
71 { "NB Elec. [W4,F]", 1 },
72 { "NB Elec. [W4-W4,V&F]", 1 },
73 { "NB Elec. [W4-W4,F]", 1 },
74 { "NB VdW & Elec. [V&F]", 1 },
75 { "NB VdW & Elec. [F]", 1 },
76 { "NB VdW & Elec. [W3,V&F]", 1 },
77 { "NB VdW & Elec. [W3,F]", 1 },
78 { "NB VdW & Elec. [W3-W3,V&F]", 1 },
79 { "NB VdW & Elec. [W3-W3,F]", 1 },
80 { "NB VdW & Elec. [W4,V&F]", 1 },
81 { "NB VdW & Elec. [W4,F]", 1 },
82 { "NB VdW & Elec. [W4-W4,V&F]", 1 },
83 { "NB VdW & Elec. [W4-W4,F]", 1 },
85 { "NB Generic kernel", 1 },
86 { "NB Generic charge grp 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 GPUs.
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 GPUs, 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. + LJ [F]", 38 }, /* nbnxn kernel LJ+RF, no ener */
103 { "NxN RF Elec. + LJ [V&F]", 54 },
104 { "NxN QSTab Elec. + LJ [F]", 41 }, /* nbnxn kernel LJ+tab, no en */
105 { "NxN QSTab Elec. + LJ [V&F]", 59 },
106 { "NxN Ewald Elec. + LJ [F]", 66 }, /* nbnxn kernel LJ+Ewald, no en */
107 { "NxN Ewald Elec. + LJ [V&F]", 107 },
108 { "NxN LJ [F]", 33 }, /* nbnxn kernel LJ, no ener */
109 { "NxN LJ [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 /* The switch function flops should be added to the LJ kernels above */
117 { "NxN LJ add F-switch [F]", 12 }, /* extra cost for LJ F-switch */
118 { "NxN LJ add F-switch [V&F]", 22 },
119 { "NxN LJ add P-switch [F]", 27 }, /* extra cost for LJ P-switch */
120 { "NxN LJ add P-switch [V&F]", 20 },
121 { "NxN LJ add LJ Ewald [F]", 36 }, /* extra cost for LJ Ewald */
122 { "NxN LJ add LJ Ewald [V&F]", 33 },
123 { "1,4 nonbonded interactions", 90 },
124 { "Born radii (Still)", 47 },
125 { "Born radii (HCT/OBC)", 183 },
126 { "Born force chain rule", 15 },
127 { "All-vs-All Still radii", 1 },
128 { "All-vs-All HCT/OBC radii", 1 },
129 { "All-vs-All Born chain rule", 1 },
130 { "Calc Weights", 36 },
132 { "Spread Q Bspline", 2 },
134 { "Gather F Bspline", 6 },
136 { "Convolution", 4 },
139 { "Reset In Box", 3 },
145 { "FENE Bonds", 58 },
146 { "Tab. Bonds", 62 },
147 { "Restraint Potential", 86 },
148 { "Linear Angles", 57 },
150 { "G96Angles", 150 },
151 { "Quartic Angles", 160 },
152 { "Tab. Angles", 169 },
154 { "Impropers", 208 },
155 { "RB-Dihedrals", 247 },
156 { "Four. Dihedrals", 247 },
157 { "Tab. Dihedrals", 227 },
158 { "Dist. Restr.", 200 },
159 { "Orient. Restr.", 200 },
160 { "Dihedral Restr.", 200 },
161 { "Pos. Restr.", 50 },
162 { "Flat-bottom posres", 50 },
163 { "Angle Restr.", 191 },
164 { "Angle Restr. Z", 164 },
165 { "Morse Potent.", 83 },
166 { "Cubic Bonds", 54 },
168 { "Polarization", 59 },
169 { "Anharmonic Polarization", 72 },
170 { "Water Pol.", 62 },
171 { "Thole Pol.", 296 },
174 { "Ext.ens. Update", 54 },
181 { "Constraint-V", 8 },
182 { "Shake-Init", 10 },
183 { "Constraint-Vir", 24 },
185 { "Virtual Site 2", 23 },
186 { "Virtual Site 3", 37 },
187 { "Virtual Site 3fd", 95 },
188 { "Virtual Site 3fad", 176 },
189 { "Virtual Site 3out", 87 },
190 { "Virtual Site 4fd", 110 },
191 { "Virtual Site 4fdn", 254 },
192 { "Virtual Site N", 15 },
193 { "Mixed Generalized Born stuff", 10 }
196 static void pr_two(FILE *out
, int c
, int i
)
200 fprintf(out
, "%c0%1d", c
, i
);
204 fprintf(out
, "%c%2d", c
, i
);
208 static void pr_difftime(FILE *out
, double dt
)
210 int ndays
, nhours
, nmins
, nsecs
;
211 gmx_bool bPrint
, bPrinted
;
213 ndays
= static_cast<int>(dt
/(24*3600));
214 dt
= dt
-24*3600*ndays
;
215 nhours
= static_cast<int>(dt
/3600);
217 nmins
= static_cast<int>(dt
/60);
219 nsecs
= static_cast<int>(dt
);
220 bPrint
= (ndays
> 0);
224 fprintf(out
, "%d", ndays
);
226 bPrint
= bPrint
|| (nhours
> 0);
231 pr_two(out
, 'd', nhours
);
235 fprintf(out
, "%d", nhours
);
238 bPrinted
= bPrinted
|| bPrint
;
239 bPrint
= bPrint
|| (nmins
> 0);
244 pr_two(out
, 'h', nmins
);
248 fprintf(out
, "%d", nmins
);
251 bPrinted
= bPrinted
|| bPrint
;
254 pr_two(out
, ':', nsecs
);
258 fprintf(out
, "%ds", nsecs
);
263 void init_nrnb(t_nrnb
*nrnb
)
267 for (i
= 0; (i
< eNRNB
); i
++)
273 void cp_nrnb(t_nrnb
*dest
, t_nrnb
*src
)
277 for (i
= 0; (i
< eNRNB
); i
++)
279 dest
->n
[i
] = src
->n
[i
];
283 void add_nrnb(t_nrnb
*dest
, t_nrnb
*s1
, t_nrnb
*s2
)
287 for (i
= 0; (i
< eNRNB
); i
++)
289 dest
->n
[i
] = s1
->n
[i
]+s2
->n
[i
];
293 void print_nrnb(FILE *out
, t_nrnb
*nrnb
)
297 for (i
= 0; (i
< eNRNB
); i
++)
301 fprintf(out
, " %-26s %10.0f.\n", nbdata
[i
].name
, nrnb
->n
[i
]);
306 void _inc_nrnb(t_nrnb
*nrnb
, int enr
, int inc
, char gmx_unused
*file
, int gmx_unused line
)
310 printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
311 nbdata
[enr
].name
, enr
, inc
, file
, line
);
315 /* Returns in enr is the index of a full nbnxn VdW kernel */
316 static gmx_bool
nrnb_is_nbnxn_vdw_kernel(int enr
)
318 return (enr
>= eNR_NBNXN_LJ_RF
&& enr
<= eNR_NBNXN_LJ_E
);
321 /* Returns in enr is the index of an nbnxn kernel addition (LJ modification) */
322 static gmx_bool
nrnb_is_nbnxn_kernel_addition(int enr
)
324 return (enr
>= eNR_NBNXN_ADD_LJ_FSW
&& enr
<= eNR_NBNXN_ADD_LJ_EWALD_E
);
327 void print_flop(FILE *out
, t_nrnb
*nrnb
, double *nbfs
, double *mflop
)
330 double mni
, frac
, tfrac
, tflop
;
331 const char *myline
= "-----------------------------------------------------------------------------";
334 for (i
= 0; (i
< eNR_NBKERNEL_ALLVSALLGB
); i
++)
336 if (std::strstr(nbdata
[i
].name
, "W3-W3") != NULL
)
338 *nbfs
+= 9e-6*nrnb
->n
[i
];
340 else if (std::strstr(nbdata
[i
].name
, "W3") != NULL
)
342 *nbfs
+= 3e-6*nrnb
->n
[i
];
344 else if (std::strstr(nbdata
[i
].name
, "W4-W4") != NULL
)
346 *nbfs
+= 10e-6*nrnb
->n
[i
];
348 else if (std::strstr(nbdata
[i
].name
, "W4") != NULL
)
350 *nbfs
+= 4e-6*nrnb
->n
[i
];
354 *nbfs
+= 1e-6*nrnb
->n
[i
];
358 for (i
= 0; (i
< eNRNB
); i
++)
360 tflop
+= 1e-6*nrnb
->n
[i
]*nbdata
[i
].flop
;
365 fprintf(out
, "No MEGA Flopsen this time\n");
370 fprintf(out
, "\n\tM E G A - F L O P S A C C O U N T I N G\n\n");
375 fprintf(out
, " NB=Group-cutoff nonbonded kernels NxN=N-by-N cluster Verlet kernels\n");
376 fprintf(out
, " RF=Reaction-Field VdW=Van der Waals QSTab=quadratic-spline table\n");
377 fprintf(out
, " W3=SPC/TIP3p W4=TIP4p (single or pairs)\n");
378 fprintf(out
, " V&F=Potential and force V=Potential only F=Force only\n\n");
380 fprintf(out
, " %-32s %16s %15s %7s\n",
381 "Computing:", "M-Number", "M-Flops", "% Flops");
382 fprintf(out
, "%s\n", myline
);
386 for (i
= 0; (i
< eNRNB
); i
++)
388 mni
= 1e-6*nrnb
->n
[i
];
389 /* Skip empty entries and nbnxn additional flops,
390 * which have been added to the kernel entry.
392 if (mni
> 0 && !nrnb_is_nbnxn_kernel_addition(i
))
396 flop
= nbdata
[i
].flop
;
397 if (nrnb_is_nbnxn_vdw_kernel(i
))
399 /* Possibly add the cost of an LJ switch/Ewald function */
400 for (j
= eNR_NBNXN_ADD_LJ_FSW
; j
<= eNR_NBNXN_ADD_LJ_EWALD
; j
+= 2)
404 /* Select the force or energy flop count */
405 e_kernel_add
= j
+ ((i
- eNR_NBNXN_LJ_RF
) % 2);
407 if (nrnb
->n
[e_kernel_add
] > 0)
409 flop
+= nbdata
[e_kernel_add
].flop
;
414 frac
= 100.0*mni
*flop
/tflop
;
418 fprintf(out
, " %-32s %16.6f %15.3f %6.1f\n",
419 nbdata
[i
].name
, mni
, mni
*flop
, frac
);
425 fprintf(out
, "%s\n", myline
);
426 fprintf(out
, " %-32s %16s %15.3f %6.1f\n",
427 "Total", "", *mflop
, tfrac
);
428 fprintf(out
, "%s\n\n", myline
);
430 if (nrnb
->n
[eNR_NBKERNEL_GENERIC
] > 0)
433 "WARNING: Using the slow generic C kernel. This is fine if you are\n"
434 "comparing different implementations or MD software. Routine\n"
435 "simulations should use a different non-bonded setup for much better\n"
441 void print_perf(FILE *out
, double time_per_thread
, double time_per_node
,
442 gmx_int64_t nsteps
, double delta_t
,
443 double nbfs
, double mflop
)
445 double wallclocktime
;
449 if (time_per_node
> 0)
451 fprintf(out
, "%12s %12s %12s %10s\n", "", "Core t (s)", "Wall t (s)", "(%)");
452 fprintf(out
, "%12s %12.3f %12.3f %10.1f\n", "Time:",
453 time_per_thread
, time_per_node
, 100.0*time_per_thread
/time_per_node
);
454 /* only print day-hour-sec format if time_per_node is more than 30 min */
455 if (time_per_node
> 30*60)
457 fprintf(out
, "%12s %12s", "", "");
458 pr_difftime(out
, time_per_node
);
462 mflop
= mflop
/time_per_node
;
463 wallclocktime
= nsteps
*delta_t
;
465 if (getenv("GMX_DETAILED_PERF_STATS") == NULL
)
467 fprintf(out
, "%12s %12s %12s\n",
468 "", "(ns/day)", "(hour/ns)");
469 fprintf(out
, "%12s %12.3f %12.3f\n", "Performance:",
470 wallclocktime
*24*3.6/time_per_node
, 1000*time_per_node
/(3600*wallclocktime
));
474 fprintf(out
, "%12s %12s %12s %12s %12s\n",
475 "", "(Mnbf/s)", (mflop
> 1000) ? "(GFlops)" : "(MFlops)",
476 "(ns/day)", "(hour/ns)");
477 fprintf(out
, "%12s %12.3f %12.3f %12.3f %12.3f\n", "Performance:",
478 nbfs
/time_per_node
, (mflop
> 1000) ? (mflop
/1000) : mflop
,
479 wallclocktime
*24*3.6/time_per_node
, 1000*time_per_node
/(3600*wallclocktime
));
484 if (getenv("GMX_DETAILED_PERF_STATS") == NULL
)
486 fprintf(out
, "%12s %14s\n",
488 fprintf(out
, "%12s %14.1f\n", "Performance:",
489 nsteps
*3600.0/time_per_node
);
493 fprintf(out
, "%12s %12s %12s %14s\n",
494 "", "(Mnbf/s)", (mflop
> 1000) ? "(GFlops)" : "(MFlops)",
496 fprintf(out
, "%12s %12.3f %12.3f %14.1f\n", "Performance:",
497 nbfs
/time_per_node
, (mflop
> 1000) ? (mflop
/1000) : mflop
,
498 nsteps
*3600.0/time_per_node
);
504 int cost_nrnb(int enr
)
506 return nbdata
[enr
].flop
;
509 const char *nrnb_str(int enr
)
511 return nbdata
[enr
].name
;
514 static const int force_index
[] = {
515 eNR_BONDS
, eNR_ANGLES
, eNR_PROPER
, eNR_IMPROPER
,
516 eNR_RB
, eNR_DISRES
, eNR_ORIRES
, eNR_POSRES
,
517 eNR_FBPOSRES
, eNR_NS
,
519 #define NFORCE_INDEX asize(force_index)
521 static const int constr_index
[] = {
522 eNR_SHAKE
, eNR_SHAKE_RIJ
, eNR_SETTLE
, eNR_UPDATE
, eNR_PCOUPL
,
523 eNR_CONSTR_VIR
, eNR_CONSTR_V
525 #define NCONSTR_INDEX asize(constr_index)
527 static double pr_av(FILE *log
, t_commrec
*cr
,
528 double fav
, double ftot
[], const char *title
)
536 fav
/= cr
->nnodes
- cr
->npmenodes
;
537 fprintf(log
, "\n %-26s", title
);
538 for (i
= 0; (i
< cr
->nnodes
); i
++)
540 dperc
= (100.0*ftot
[i
])/fav
;
541 unb
= std::max(unb
, dperc
);
542 perc
= static_cast<int>(dperc
);
543 fprintf(log
, "%3d ", perc
);
547 perc
= static_cast<int>(10000.0/unb
);
548 fprintf(log
, "%6d%%\n\n", perc
);
552 fprintf(log
, "\n\n");
558 void pr_load(FILE *log
, t_commrec
*cr
, t_nrnb nrnb
[])
561 double dperc
, unb
, uf
, us
;
567 snew(ftot
, cr
->nnodes
);
568 snew(stot
, cr
->nnodes
);
570 for (i
= 0; (i
< cr
->nnodes
); i
++)
572 add_nrnb(av
, av
, &(nrnb
[i
]));
573 /* Cost due to forces */
574 for (j
= 0; (j
< eNR_NBKERNEL_ALLVSALLGB
); j
++)
576 ftot
[i
] += nrnb
[i
].n
[j
]*cost_nrnb(j
);
578 for (j
= 0; (j
< NFORCE_INDEX
); j
++)
580 ftot
[i
] += nrnb
[i
].n
[force_index
[j
]]*cost_nrnb(force_index
[j
]);
583 for (j
= 0; (j
< NCONSTR_INDEX
); j
++)
585 stot
[i
] += nrnb
[i
].n
[constr_index
[j
]]*cost_nrnb(constr_index
[j
]);
588 for (j
= 0; (j
< eNRNB
); j
++)
590 av
->n
[j
] = av
->n
[j
]/static_cast<double>(cr
->nnodes
- cr
->npmenodes
);
593 fprintf(log
, "\nDetailed load balancing info in percentage of average\n");
595 fprintf(log
, " Type RANK:");
596 for (i
= 0; (i
< cr
->nnodes
); i
++)
598 fprintf(log
, "%3d ", i
);
600 fprintf(log
, "Scaling\n");
601 fprintf(log
, "---------------------------");
602 for (i
= 0; (i
< cr
->nnodes
); i
++)
604 fprintf(log
, "----");
606 fprintf(log
, "-------\n");
608 for (j
= 0; (j
< eNRNB
); j
++)
613 fprintf(log
, " %-26s", nrnb_str(j
));
614 for (i
= 0; (i
< cr
->nnodes
); i
++)
616 dperc
= (100.0*nrnb
[i
].n
[j
])/av
->n
[j
];
617 unb
= std::max(unb
, dperc
);
618 perc
= static_cast<int>(dperc
);
619 fprintf(log
, "%3d ", perc
);
623 perc
= static_cast<int>(10000.0/unb
);
624 fprintf(log
, "%6d%%\n", perc
);
633 for (i
= 0; (i
< cr
->nnodes
); i
++)
638 uf
= pr_av(log
, cr
, fav
, ftot
, "Total Force");
639 us
= pr_av(log
, cr
, sav
, stot
, "Total Constr.");
641 unb
= (uf
*fav
+us
*sav
)/(fav
+sav
);
645 fprintf(log
, "\nTotal Scaling: %.0f%% of max performance\n\n", unb
);