Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / gmxlib / nrnb.cpp
blob08de8f3af1e8d09f66c79559340f31f7e9fcd81d
1 /*
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,2019, 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.
37 #include "gmxpre.h"
39 #include "nrnb.h"
41 #include <cstdlib>
42 #include <cstring>
44 #include <algorithm>
45 #include <vector>
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/md_enums.h"
49 #include "gromacs/utility/arraysize.h"
51 typedef struct {
52 const char *name;
53 int flop;
54 } t_nrnb_data;
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 },
62 { "NB VdW [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 Free energy kernel", 1 },
88 { "Pair Search distance check", 9 }, /* nbnxn pair dist. check */
89 /* nbnxn kernel flops are based on inner-loops without exclusion checks.
90 * Plain Coulomb runs through the RF kernels, except with GPUs.
91 * invsqrt is counted as 6 flops: 1 for _mm_rsqt_ps + 5 for iteration.
92 * The flops are equal for plain-C, x86 SIMD and GPUs, except for:
93 * - plain-C kernel uses one flop more for Coulomb-only (F) than listed
94 * - x86 SIMD LJ geom-comb.rule kernels (fastest) use 2 more flops
95 * - x86 SIMD LJ LB-comb.rule kernels (fast) use 3 (8 for F+E) more flops
96 * - GPU always does exclusions, which requires 2-4 flops, but as invsqrt
97 * is always counted as 6 flops, this roughly compensates.
99 { "NxN RF Elec. + LJ [F]", 38 }, /* nbnxn kernel LJ+RF, no ener */
100 { "NxN RF Elec. + LJ [V&F]", 54 },
101 { "NxN QSTab Elec. + LJ [F]", 41 }, /* nbnxn kernel LJ+tab, no en */
102 { "NxN QSTab Elec. + LJ [V&F]", 59 },
103 { "NxN Ewald Elec. + LJ [F]", 66 }, /* nbnxn kernel LJ+Ewald, no en */
104 { "NxN Ewald Elec. + LJ [V&F]", 107 },
105 { "NxN LJ [F]", 33 }, /* nbnxn kernel LJ, no ener */
106 { "NxN LJ [V&F]", 43 },
107 { "NxN RF Electrostatics [F]", 31 }, /* nbnxn kernel RF, no ener */
108 { "NxN RF Electrostatics [V&F]", 36 },
109 { "NxN QSTab Elec. [F]", 34 }, /* nbnxn kernel tab, no ener */
110 { "NxN QSTab Elec. [V&F]", 41 },
111 { "NxN Ewald Elec. [F]", 61 }, /* nbnxn kernel Ewald, no ener */
112 { "NxN Ewald Elec. [V&F]", 84 },
113 /* The switch function flops should be added to the LJ kernels above */
114 { "NxN LJ add F-switch [F]", 12 }, /* extra cost for LJ F-switch */
115 { "NxN LJ add F-switch [V&F]", 22 },
116 { "NxN LJ add P-switch [F]", 27 }, /* extra cost for LJ P-switch */
117 { "NxN LJ add P-switch [V&F]", 20 },
118 { "NxN LJ add LJ Ewald [F]", 36 }, /* extra cost for LJ Ewald */
119 { "NxN LJ add LJ Ewald [V&F]", 33 },
120 { "1,4 nonbonded interactions", 90 },
121 { "Calc Weights", 36 },
122 { "Spread Q", 6 },
123 { "Spread Q Bspline", 2 },
124 { "Gather F", 23 },
125 { "Gather F Bspline", 6 },
126 { "3D-FFT", 8 },
127 { "Convolution", 4 },
128 { "Solve PME", 64 },
129 { "NS-Pairs", 21 },
130 { "Reset In Box", 3 },
131 { "Shift-X", 6 },
132 { "CG-CoM", 3 },
133 { "Sum Forces", 1 },
134 { "Bonds", 59 },
135 { "G96Bonds", 44 },
136 { "FENE Bonds", 58 },
137 { "Tab. Bonds", 62 },
138 { "Restraint Potential", 86 },
139 { "Linear Angles", 57 },
140 { "Angles", 168 },
141 { "G96Angles", 150 },
142 { "Quartic Angles", 160 },
143 { "Tab. Angles", 169 },
144 { "Propers", 229 },
145 { "Impropers", 208 },
146 { "RB-Dihedrals", 247 },
147 { "Four. Dihedrals", 247 },
148 { "Tab. Dihedrals", 227 },
149 { "Dist. Restr.", 200 },
150 { "Orient. Restr.", 200 },
151 { "Dihedral Restr.", 200 },
152 { "Pos. Restr.", 50 },
153 { "Flat-bottom posres", 50 },
154 { "Angle Restr.", 191 },
155 { "Angle Restr. Z", 164 },
156 { "Morse Potent.", 83 },
157 { "Cubic Bonds", 54 },
158 { "Walls", 31 },
159 { "Polarization", 59 },
160 { "Anharmonic Polarization", 72 },
161 { "Water Pol.", 62 },
162 { "Thole Pol.", 296 },
163 { "Virial", 18 },
164 { "Update", 31 },
165 { "Ext.ens. Update", 54 },
166 { "Stop-CM", 10 },
167 { "P-Coupling", 6 },
168 { "Calc-Ekin", 27 },
169 { "Lincs", 60 },
170 { "Lincs-Mat", 4 },
171 { "Shake", 30 },
172 { "Constraint-V", 8 },
173 { "Shake-Init", 10 },
174 { "Constraint-Vir", 24 },
175 { "Settle", 323 },
176 { "Virtual Site 2", 23 },
177 { "Virtual Site 3", 37 },
178 { "Virtual Site 3fd", 95 },
179 { "Virtual Site 3fad", 176 },
180 { "Virtual Site 3out", 87 },
181 { "Virtual Site 4fd", 110 },
182 { "Virtual Site 4fdn", 254 },
183 { "Virtual Site N", 15 },
184 { "CMAP", 1700 }, // Estimate!
185 { "Urey-Bradley", 183 },
186 { "Cross-Bond-Bond", 163 },
187 { "Cross-Bond-Angle", 163 }
190 static void pr_two(FILE *out, int c, int i)
192 if (i < 10)
194 fprintf(out, "%c0%1d", c, i);
196 else
198 fprintf(out, "%c%2d", c, i);
202 static void pr_difftime(FILE *out, double dt)
204 int ndays, nhours, nmins, nsecs;
205 gmx_bool bPrint, bPrinted;
207 ndays = static_cast<int>(dt/(24*3600));
208 dt = dt-24*3600*ndays;
209 nhours = static_cast<int>(dt/3600);
210 dt = dt-3600*nhours;
211 nmins = static_cast<int>(dt/60);
212 dt = dt-nmins*60;
213 nsecs = static_cast<int>(dt);
214 bPrint = (ndays > 0);
215 bPrinted = bPrint;
216 if (bPrint)
218 fprintf(out, "%d", ndays);
220 bPrint = bPrint || (nhours > 0);
221 if (bPrint)
223 if (bPrinted)
225 pr_two(out, 'd', nhours);
227 else
229 fprintf(out, "%d", nhours);
232 bPrinted = bPrinted || bPrint;
233 bPrint = bPrint || (nmins > 0);
234 if (bPrint)
236 if (bPrinted)
238 pr_two(out, 'h', nmins);
240 else
242 fprintf(out, "%d", nmins);
245 bPrinted = bPrinted || bPrint;
246 if (bPrinted)
248 pr_two(out, ':', nsecs);
250 else
252 fprintf(out, "%ds", nsecs);
254 fprintf(out, "\n");
257 void clear_nrnb(t_nrnb *nrnb)
259 int i;
261 for (i = 0; (i < eNRNB); i++)
263 nrnb->n[i] = 0.0;
267 void add_nrnb(t_nrnb *dest, t_nrnb *s1, t_nrnb *s2)
269 int i;
271 for (i = 0; (i < eNRNB); i++)
273 dest->n[i] = s1->n[i]+s2->n[i];
277 void print_nrnb(FILE *out, t_nrnb *nrnb)
279 int i;
281 for (i = 0; (i < eNRNB); i++)
283 if (nrnb->n[i] > 0)
285 fprintf(out, " %-26s %10.0f.\n", nbdata[i].name, nrnb->n[i]);
290 void _inc_nrnb(t_nrnb *nrnb, int enr, int inc, char gmx_unused *file, int gmx_unused line)
292 nrnb->n[enr] += inc;
293 #ifdef DEBUG_NRNB
294 printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
295 nbdata[enr].name, enr, inc, file, line);
296 #endif
299 /* Returns in enr is the index of a full nbnxn VdW kernel */
300 static gmx_bool nrnb_is_nbnxn_vdw_kernel(int enr)
302 return (enr >= eNR_NBNXN_LJ_RF && enr <= eNR_NBNXN_LJ_E);
305 /* Returns in enr is the index of an nbnxn kernel addition (LJ modification) */
306 static gmx_bool nrnb_is_nbnxn_kernel_addition(int enr)
308 return (enr >= eNR_NBNXN_ADD_LJ_FSW && enr <= eNR_NBNXN_ADD_LJ_EWALD_E);
311 void print_flop(FILE *out, t_nrnb *nrnb, double *nbfs, double *mflop)
313 int i, j;
314 double mni, frac, tfrac, tflop;
315 const char *myline = "-----------------------------------------------------------------------------";
317 *nbfs = 0.0;
318 for (i = 0; (i < eNR_NBKERNEL_TOTAL_NR); i++)
320 if (std::strstr(nbdata[i].name, "W3-W3") != nullptr)
322 *nbfs += 9e-6*nrnb->n[i];
324 else if (std::strstr(nbdata[i].name, "W3") != nullptr)
326 *nbfs += 3e-6*nrnb->n[i];
328 else if (std::strstr(nbdata[i].name, "W4-W4") != nullptr)
330 *nbfs += 10e-6*nrnb->n[i];
332 else if (std::strstr(nbdata[i].name, "W4") != nullptr)
334 *nbfs += 4e-6*nrnb->n[i];
336 else
338 *nbfs += 1e-6*nrnb->n[i];
341 tflop = 0;
342 for (i = 0; (i < eNRNB); i++)
344 tflop += 1e-6*nrnb->n[i]*nbdata[i].flop;
347 if (tflop == 0)
349 fprintf(out, "No MEGA Flopsen this time\n");
350 return;
352 if (out)
354 fprintf(out, "\n\tM E G A - F L O P S A C C O U N T I N G\n\n");
357 if (out)
359 fprintf(out, " NB=Group-cutoff nonbonded kernels NxN=N-by-N cluster Verlet kernels\n");
360 fprintf(out, " RF=Reaction-Field VdW=Van der Waals QSTab=quadratic-spline table\n");
361 fprintf(out, " W3=SPC/TIP3p W4=TIP4p (single or pairs)\n");
362 fprintf(out, " V&F=Potential and force V=Potential only F=Force only\n\n");
364 fprintf(out, " %-32s %16s %15s %7s\n",
365 "Computing:", "M-Number", "M-Flops", "% Flops");
366 fprintf(out, "%s\n", myline);
368 *mflop = 0.0;
369 tfrac = 0.0;
370 for (i = 0; (i < eNRNB); i++)
372 mni = 1e-6*nrnb->n[i];
373 /* Skip empty entries and nbnxn additional flops,
374 * which have been added to the kernel entry.
376 if (mni > 0 && !nrnb_is_nbnxn_kernel_addition(i))
378 int flop;
380 flop = nbdata[i].flop;
381 if (nrnb_is_nbnxn_vdw_kernel(i))
383 /* Possibly add the cost of an LJ switch/Ewald function */
384 for (j = eNR_NBNXN_ADD_LJ_FSW; j <= eNR_NBNXN_ADD_LJ_EWALD; j += 2)
386 int e_kernel_add;
388 /* Select the force or energy flop count */
389 e_kernel_add = j + ((i - eNR_NBNXN_LJ_RF) % 2);
391 if (nrnb->n[e_kernel_add] > 0)
393 flop += nbdata[e_kernel_add].flop;
397 *mflop += mni*flop;
398 frac = 100.0*mni*flop/tflop;
399 tfrac += frac;
400 if (out != nullptr)
402 fprintf(out, " %-32s %16.6f %15.3f %6.1f\n",
403 nbdata[i].name, mni, mni*flop, frac);
407 if (out)
409 fprintf(out, "%s\n", myline);
410 fprintf(out, " %-32s %16s %15.3f %6.1f\n",
411 "Total", "", *mflop, tfrac);
412 fprintf(out, "%s\n\n", myline);
414 if (nrnb->n[eNR_NBKERNEL_GENERIC] > 0)
416 fprintf(out,
417 "WARNING: Using the slow generic C kernel. This is fine if you are\n"
418 "comparing different implementations or MD software. Routine\n"
419 "simulations should use a different non-bonded setup for much better\n"
420 "performance.\n\n");
425 void print_perf(FILE *out, double time_per_thread, double time_per_node,
426 int64_t nsteps, double delta_t,
427 double nbfs, double mflop)
429 double wallclocktime;
431 fprintf(out, "\n");
433 if (time_per_node > 0)
435 fprintf(out, "%12s %12s %12s %10s\n", "", "Core t (s)", "Wall t (s)", "(%)");
436 fprintf(out, "%12s %12.3f %12.3f %10.1f\n", "Time:",
437 time_per_thread, time_per_node, 100.0*time_per_thread/time_per_node);
438 /* only print day-hour-sec format if time_per_node is more than 30 min */
439 if (time_per_node > 30*60)
441 fprintf(out, "%12s %12s", "", "");
442 pr_difftime(out, time_per_node);
444 if (delta_t > 0)
446 mflop = mflop/time_per_node;
447 wallclocktime = nsteps*delta_t;
449 if (getenv("GMX_DETAILED_PERF_STATS") == nullptr)
451 fprintf(out, "%12s %12s %12s\n",
452 "", "(ns/day)", "(hour/ns)");
453 fprintf(out, "%12s %12.3f %12.3f\n", "Performance:",
454 wallclocktime*24*3.6/time_per_node, 1000*time_per_node/(3600*wallclocktime));
456 else
458 fprintf(out, "%12s %12s %12s %12s %12s\n",
459 "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
460 "(ns/day)", "(hour/ns)");
461 fprintf(out, "%12s %12.3f %12.3f %12.3f %12.3f\n", "Performance:",
462 nbfs/time_per_node, (mflop > 1000) ? (mflop/1000) : mflop,
463 wallclocktime*24*3.6/time_per_node, 1000*time_per_node/(3600*wallclocktime));
466 else
468 if (getenv("GMX_DETAILED_PERF_STATS") == nullptr)
470 fprintf(out, "%12s %14s\n",
471 "", "(steps/hour)");
472 fprintf(out, "%12s %14.1f\n", "Performance:",
473 nsteps*3600.0/time_per_node);
475 else
477 fprintf(out, "%12s %12s %12s %14s\n",
478 "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
479 "(steps/hour)");
480 fprintf(out, "%12s %12.3f %12.3f %14.1f\n", "Performance:",
481 nbfs/time_per_node, (mflop > 1000) ? (mflop/1000) : mflop,
482 nsteps*3600.0/time_per_node);
488 int cost_nrnb(int enr)
490 return nbdata[enr].flop;
493 const char *nrnb_str(int enr)
495 return nbdata[enr].name;
498 static const int force_index[] = {
499 eNR_BONDS, eNR_ANGLES, eNR_PROPER, eNR_IMPROPER,
500 eNR_RB, eNR_DISRES, eNR_ORIRES, eNR_POSRES,
501 eNR_FBPOSRES, eNR_NS,
503 #define NFORCE_INDEX asize(force_index)
505 static const int constr_index[] = {
506 eNR_SHAKE, eNR_SHAKE_RIJ, eNR_SETTLE, eNR_UPDATE, eNR_PCOUPL,
507 eNR_CONSTR_VIR, eNR_CONSTR_V
509 #define NCONSTR_INDEX asize(constr_index)
511 static double pr_av(FILE *log, t_commrec *cr,
512 double fav, const double ftot[], const char *title)
514 int i, perc;
515 double dperc, unb;
517 unb = 0;
518 if (fav > 0)
520 fav /= cr->nnodes - cr->npmenodes;
521 fprintf(log, "\n %-26s", title);
522 for (i = 0; (i < cr->nnodes); i++)
524 dperc = (100.0*ftot[i])/fav;
525 unb = std::max(unb, dperc);
526 perc = static_cast<int>(dperc);
527 fprintf(log, "%3d ", perc);
529 if (unb > 0)
531 perc = static_cast<int>(10000.0/unb);
532 fprintf(log, "%6d%%\n\n", perc);
534 else
536 fprintf(log, "\n\n");
539 return unb;
542 void pr_load(FILE *log, t_commrec *cr, t_nrnb nrnb[])
544 t_nrnb av;
546 std::vector<double> ftot(cr->nnodes);
547 std::vector<double> stot(cr->nnodes);
548 for (int i = 0; (i < cr->nnodes); i++)
550 add_nrnb(&av, &av, &(nrnb[i]));
551 /* Cost due to forces */
552 for (int j = 0; (j < eNR_NBKERNEL_TOTAL_NR); j++)
554 ftot[i] += nrnb[i].n[j]*cost_nrnb(j);
556 for (int j = 0; (j < NFORCE_INDEX); j++)
558 ftot[i] += nrnb[i].n[force_index[j]]*cost_nrnb(force_index[j]);
560 /* Due to shake */
561 for (int j = 0; (j < NCONSTR_INDEX); j++)
563 stot[i] += nrnb[i].n[constr_index[j]]*cost_nrnb(constr_index[j]);
566 for (int j = 0; (j < eNRNB); j++)
568 av.n[j] = av.n[j]/static_cast<double>(cr->nnodes - cr->npmenodes);
571 fprintf(log, "\nDetailed load balancing info in percentage of average\n");
573 fprintf(log, " Type RANK:");
574 for (int i = 0; (i < cr->nnodes); i++)
576 fprintf(log, "%3d ", i);
578 fprintf(log, "Scaling\n");
579 fprintf(log, "---------------------------");
580 for (int i = 0; (i < cr->nnodes); i++)
582 fprintf(log, "----");
584 fprintf(log, "-------\n");
586 for (int j = 0; (j < eNRNB); j++)
588 double unb = 100.0;
589 if (av.n[j] > 0)
591 double dperc;
592 int perc;
593 fprintf(log, " %-26s", nrnb_str(j));
594 for (int i = 0; (i < cr->nnodes); i++)
596 dperc = (100.0*nrnb[i].n[j])/av.n[j];
597 unb = std::max(unb, dperc);
598 perc = static_cast<int>(dperc);
599 fprintf(log, "%3d ", perc);
601 if (unb > 0)
603 perc = static_cast<int>(10000.0/unb);
604 fprintf(log, "%6d%%\n", perc);
606 else
608 fprintf(log, "\n");
612 double fav = 0;
613 double sav = 0;
614 for (int i = 0; (i < cr->nnodes); i++)
616 fav += ftot[i];
617 sav += stot[i];
619 double uf = pr_av(log, cr, fav, ftot.data(), "Total Force");
620 double us = pr_av(log, cr, sav, stot.data(), "Total Constr.");
622 double unb = (uf*fav+us*sav)/(fav+sav);
623 if (unb > 0)
625 unb = 10000.0/unb;
626 fprintf(log, "\nTotal Scaling: %.0f%% of max performance\n\n", unb);