prepareGpuKernelArguments() and launchGpuKernel() are added
[gromacs.git] / src / gromacs / gmxlib / nrnb.cpp
blobba5ec1da26525a10d2abfb9ad3d5fa58018cbf6a
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, 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>
46 #include "gromacs/mdtypes/commrec.h"
47 #include "gromacs/mdtypes/md_enums.h"
48 #include "gromacs/utility/arraysize.h"
49 #include "gromacs/utility/smalloc.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 },
87 { "NB All-vs-all", 1 },
89 { "Pair Search distance check", 9 }, /* nbnxn pair dist. check */
90 /* nbnxn kernel flops are based on inner-loops without exclusion checks.
91 * Plain Coulomb runs through the RF kernels, except with GPUs.
92 * invsqrt is counted as 6 flops: 1 for _mm_rsqt_ps + 5 for iteration.
93 * The flops are equal for plain-C, x86 SIMD and GPUs, except for:
94 * - plain-C kernel uses one flop more for Coulomb-only (F) than listed
95 * - x86 SIMD LJ geom-comb.rule kernels (fastest) use 2 more flops
96 * - x86 SIMD LJ LB-comb.rule kernels (fast) use 3 (8 for F+E) more flops
97 * - GPU always does exclusions, which requires 2-4 flops, but as invsqrt
98 * is always counted as 6 flops, this roughly compensates.
100 { "NxN RF Elec. + LJ [F]", 38 }, /* nbnxn kernel LJ+RF, no ener */
101 { "NxN RF Elec. + LJ [V&F]", 54 },
102 { "NxN QSTab Elec. + LJ [F]", 41 }, /* nbnxn kernel LJ+tab, no en */
103 { "NxN QSTab Elec. + LJ [V&F]", 59 },
104 { "NxN Ewald Elec. + LJ [F]", 66 }, /* nbnxn kernel LJ+Ewald, no en */
105 { "NxN Ewald Elec. + LJ [V&F]", 107 },
106 { "NxN LJ [F]", 33 }, /* nbnxn kernel LJ, no ener */
107 { "NxN LJ [V&F]", 43 },
108 { "NxN RF Electrostatics [F]", 31 }, /* nbnxn kernel RF, no ener */
109 { "NxN RF Electrostatics [V&F]", 36 },
110 { "NxN QSTab Elec. [F]", 34 }, /* nbnxn kernel tab, no ener */
111 { "NxN QSTab Elec. [V&F]", 41 },
112 { "NxN Ewald Elec. [F]", 61 }, /* nbnxn kernel Ewald, no ener */
113 { "NxN Ewald Elec. [V&F]", 84 },
114 /* The switch function flops should be added to the LJ kernels above */
115 { "NxN LJ add F-switch [F]", 12 }, /* extra cost for LJ F-switch */
116 { "NxN LJ add F-switch [V&F]", 22 },
117 { "NxN LJ add P-switch [F]", 27 }, /* extra cost for LJ P-switch */
118 { "NxN LJ add P-switch [V&F]", 20 },
119 { "NxN LJ add LJ Ewald [F]", 36 }, /* extra cost for LJ Ewald */
120 { "NxN LJ add LJ Ewald [V&F]", 33 },
121 { "1,4 nonbonded interactions", 90 },
122 { "Calc Weights", 36 },
123 { "Spread Q", 6 },
124 { "Spread Q Bspline", 2 },
125 { "Gather F", 23 },
126 { "Gather F Bspline", 6 },
127 { "3D-FFT", 8 },
128 { "Convolution", 4 },
129 { "Solve PME", 64 },
130 { "NS-Pairs", 21 },
131 { "Reset In Box", 3 },
132 { "Shift-X", 6 },
133 { "CG-CoM", 3 },
134 { "Sum Forces", 1 },
135 { "Bonds", 59 },
136 { "G96Bonds", 44 },
137 { "FENE Bonds", 58 },
138 { "Tab. Bonds", 62 },
139 { "Restraint Potential", 86 },
140 { "Linear Angles", 57 },
141 { "Angles", 168 },
142 { "G96Angles", 150 },
143 { "Quartic Angles", 160 },
144 { "Tab. Angles", 169 },
145 { "Propers", 229 },
146 { "Impropers", 208 },
147 { "RB-Dihedrals", 247 },
148 { "Four. Dihedrals", 247 },
149 { "Tab. Dihedrals", 227 },
150 { "Dist. Restr.", 200 },
151 { "Orient. Restr.", 200 },
152 { "Dihedral Restr.", 200 },
153 { "Pos. Restr.", 50 },
154 { "Flat-bottom posres", 50 },
155 { "Angle Restr.", 191 },
156 { "Angle Restr. Z", 164 },
157 { "Morse Potent.", 83 },
158 { "Cubic Bonds", 54 },
159 { "Walls", 31 },
160 { "Polarization", 59 },
161 { "Anharmonic Polarization", 72 },
162 { "Water Pol.", 62 },
163 { "Thole Pol.", 296 },
164 { "Virial", 18 },
165 { "Update", 31 },
166 { "Ext.ens. Update", 54 },
167 { "Stop-CM", 10 },
168 { "P-Coupling", 6 },
169 { "Calc-Ekin", 27 },
170 { "Lincs", 60 },
171 { "Lincs-Mat", 4 },
172 { "Shake", 30 },
173 { "Constraint-V", 8 },
174 { "Shake-Init", 10 },
175 { "Constraint-Vir", 24 },
176 { "Settle", 323 },
177 { "Virtual Site 2", 23 },
178 { "Virtual Site 3", 37 },
179 { "Virtual Site 3fd", 95 },
180 { "Virtual Site 3fad", 176 },
181 { "Virtual Site 3out", 87 },
182 { "Virtual Site 4fd", 110 },
183 { "Virtual Site 4fdn", 254 },
184 { "Virtual Site N", 15 },
185 { "CMAP", 1700 }, // Estimate!
186 { "Urey-Bradley", 183 },
187 { "Cross-Bond-Bond", 163 },
188 { "Cross-Bond-Angle", 163 }
191 static void pr_two(FILE *out, int c, int i)
193 if (i < 10)
195 fprintf(out, "%c0%1d", c, i);
197 else
199 fprintf(out, "%c%2d", c, i);
203 static void pr_difftime(FILE *out, double dt)
205 int ndays, nhours, nmins, nsecs;
206 gmx_bool bPrint, bPrinted;
208 ndays = static_cast<int>(dt/(24*3600));
209 dt = dt-24*3600*ndays;
210 nhours = static_cast<int>(dt/3600);
211 dt = dt-3600*nhours;
212 nmins = static_cast<int>(dt/60);
213 dt = dt-nmins*60;
214 nsecs = static_cast<int>(dt);
215 bPrint = (ndays > 0);
216 bPrinted = bPrint;
217 if (bPrint)
219 fprintf(out, "%d", ndays);
221 bPrint = bPrint || (nhours > 0);
222 if (bPrint)
224 if (bPrinted)
226 pr_two(out, 'd', nhours);
228 else
230 fprintf(out, "%d", nhours);
233 bPrinted = bPrinted || bPrint;
234 bPrint = bPrint || (nmins > 0);
235 if (bPrint)
237 if (bPrinted)
239 pr_two(out, 'h', nmins);
241 else
243 fprintf(out, "%d", nmins);
246 bPrinted = bPrinted || bPrint;
247 if (bPrinted)
249 pr_two(out, ':', nsecs);
251 else
253 fprintf(out, "%ds", nsecs);
255 fprintf(out, "\n");
258 void init_nrnb(t_nrnb *nrnb)
260 int i;
262 for (i = 0; (i < eNRNB); i++)
264 nrnb->n[i] = 0.0;
268 void cp_nrnb(t_nrnb *dest, t_nrnb *src)
270 int i;
272 for (i = 0; (i < eNRNB); i++)
274 dest->n[i] = src->n[i];
278 void add_nrnb(t_nrnb *dest, t_nrnb *s1, t_nrnb *s2)
280 int i;
282 for (i = 0; (i < eNRNB); i++)
284 dest->n[i] = s1->n[i]+s2->n[i];
288 void print_nrnb(FILE *out, t_nrnb *nrnb)
290 int i;
292 for (i = 0; (i < eNRNB); i++)
294 if (nrnb->n[i] > 0)
296 fprintf(out, " %-26s %10.0f.\n", nbdata[i].name, nrnb->n[i]);
301 void _inc_nrnb(t_nrnb *nrnb, int enr, int inc, char gmx_unused *file, int gmx_unused line)
303 nrnb->n[enr] += inc;
304 #ifdef DEBUG_NRNB
305 printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
306 nbdata[enr].name, enr, inc, file, line);
307 #endif
310 /* Returns in enr is the index of a full nbnxn VdW kernel */
311 static gmx_bool nrnb_is_nbnxn_vdw_kernel(int enr)
313 return (enr >= eNR_NBNXN_LJ_RF && enr <= eNR_NBNXN_LJ_E);
316 /* Returns in enr is the index of an nbnxn kernel addition (LJ modification) */
317 static gmx_bool nrnb_is_nbnxn_kernel_addition(int enr)
319 return (enr >= eNR_NBNXN_ADD_LJ_FSW && enr <= eNR_NBNXN_ADD_LJ_EWALD_E);
322 void print_flop(FILE *out, t_nrnb *nrnb, double *nbfs, double *mflop)
324 int i, j;
325 double mni, frac, tfrac, tflop;
326 const char *myline = "-----------------------------------------------------------------------------";
328 *nbfs = 0.0;
329 for (i = 0; (i < eNR_NBKERNEL_TOTAL_NR); i++)
331 if (std::strstr(nbdata[i].name, "W3-W3") != nullptr)
333 *nbfs += 9e-6*nrnb->n[i];
335 else if (std::strstr(nbdata[i].name, "W3") != nullptr)
337 *nbfs += 3e-6*nrnb->n[i];
339 else if (std::strstr(nbdata[i].name, "W4-W4") != nullptr)
341 *nbfs += 10e-6*nrnb->n[i];
343 else if (std::strstr(nbdata[i].name, "W4") != nullptr)
345 *nbfs += 4e-6*nrnb->n[i];
347 else
349 *nbfs += 1e-6*nrnb->n[i];
352 tflop = 0;
353 for (i = 0; (i < eNRNB); i++)
355 tflop += 1e-6*nrnb->n[i]*nbdata[i].flop;
358 if (tflop == 0)
360 fprintf(out, "No MEGA Flopsen this time\n");
361 return;
363 if (out)
365 fprintf(out, "\n\tM E G A - F L O P S A C C O U N T I N G\n\n");
368 if (out)
370 fprintf(out, " NB=Group-cutoff nonbonded kernels NxN=N-by-N cluster Verlet kernels\n");
371 fprintf(out, " RF=Reaction-Field VdW=Van der Waals QSTab=quadratic-spline table\n");
372 fprintf(out, " W3=SPC/TIP3p W4=TIP4p (single or pairs)\n");
373 fprintf(out, " V&F=Potential and force V=Potential only F=Force only\n\n");
375 fprintf(out, " %-32s %16s %15s %7s\n",
376 "Computing:", "M-Number", "M-Flops", "% Flops");
377 fprintf(out, "%s\n", myline);
379 *mflop = 0.0;
380 tfrac = 0.0;
381 for (i = 0; (i < eNRNB); i++)
383 mni = 1e-6*nrnb->n[i];
384 /* Skip empty entries and nbnxn additional flops,
385 * which have been added to the kernel entry.
387 if (mni > 0 && !nrnb_is_nbnxn_kernel_addition(i))
389 int flop;
391 flop = nbdata[i].flop;
392 if (nrnb_is_nbnxn_vdw_kernel(i))
394 /* Possibly add the cost of an LJ switch/Ewald function */
395 for (j = eNR_NBNXN_ADD_LJ_FSW; j <= eNR_NBNXN_ADD_LJ_EWALD; j += 2)
397 int e_kernel_add;
399 /* Select the force or energy flop count */
400 e_kernel_add = j + ((i - eNR_NBNXN_LJ_RF) % 2);
402 if (nrnb->n[e_kernel_add] > 0)
404 flop += nbdata[e_kernel_add].flop;
408 *mflop += mni*flop;
409 frac = 100.0*mni*flop/tflop;
410 tfrac += frac;
411 if (out != nullptr)
413 fprintf(out, " %-32s %16.6f %15.3f %6.1f\n",
414 nbdata[i].name, mni, mni*flop, frac);
418 if (out)
420 fprintf(out, "%s\n", myline);
421 fprintf(out, " %-32s %16s %15.3f %6.1f\n",
422 "Total", "", *mflop, tfrac);
423 fprintf(out, "%s\n\n", myline);
425 if (nrnb->n[eNR_NBKERNEL_GENERIC] > 0)
427 fprintf(out,
428 "WARNING: Using the slow generic C kernel. This is fine if you are\n"
429 "comparing different implementations or MD software. Routine\n"
430 "simulations should use a different non-bonded setup for much better\n"
431 "performance.\n\n");
436 void print_perf(FILE *out, double time_per_thread, double time_per_node,
437 gmx_int64_t nsteps, double delta_t,
438 double nbfs, double mflop)
440 double wallclocktime;
442 fprintf(out, "\n");
444 if (time_per_node > 0)
446 fprintf(out, "%12s %12s %12s %10s\n", "", "Core t (s)", "Wall t (s)", "(%)");
447 fprintf(out, "%12s %12.3f %12.3f %10.1f\n", "Time:",
448 time_per_thread, time_per_node, 100.0*time_per_thread/time_per_node);
449 /* only print day-hour-sec format if time_per_node is more than 30 min */
450 if (time_per_node > 30*60)
452 fprintf(out, "%12s %12s", "", "");
453 pr_difftime(out, time_per_node);
455 if (delta_t > 0)
457 mflop = mflop/time_per_node;
458 wallclocktime = nsteps*delta_t;
460 if (getenv("GMX_DETAILED_PERF_STATS") == nullptr)
462 fprintf(out, "%12s %12s %12s\n",
463 "", "(ns/day)", "(hour/ns)");
464 fprintf(out, "%12s %12.3f %12.3f\n", "Performance:",
465 wallclocktime*24*3.6/time_per_node, 1000*time_per_node/(3600*wallclocktime));
467 else
469 fprintf(out, "%12s %12s %12s %12s %12s\n",
470 "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
471 "(ns/day)", "(hour/ns)");
472 fprintf(out, "%12s %12.3f %12.3f %12.3f %12.3f\n", "Performance:",
473 nbfs/time_per_node, (mflop > 1000) ? (mflop/1000) : mflop,
474 wallclocktime*24*3.6/time_per_node, 1000*time_per_node/(3600*wallclocktime));
477 else
479 if (getenv("GMX_DETAILED_PERF_STATS") == nullptr)
481 fprintf(out, "%12s %14s\n",
482 "", "(steps/hour)");
483 fprintf(out, "%12s %14.1f\n", "Performance:",
484 nsteps*3600.0/time_per_node);
486 else
488 fprintf(out, "%12s %12s %12s %14s\n",
489 "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
490 "(steps/hour)");
491 fprintf(out, "%12s %12.3f %12.3f %14.1f\n", "Performance:",
492 nbfs/time_per_node, (mflop > 1000) ? (mflop/1000) : mflop,
493 nsteps*3600.0/time_per_node);
499 int cost_nrnb(int enr)
501 return nbdata[enr].flop;
504 const char *nrnb_str(int enr)
506 return nbdata[enr].name;
509 static const int force_index[] = {
510 eNR_BONDS, eNR_ANGLES, eNR_PROPER, eNR_IMPROPER,
511 eNR_RB, eNR_DISRES, eNR_ORIRES, eNR_POSRES,
512 eNR_FBPOSRES, eNR_NS,
514 #define NFORCE_INDEX asize(force_index)
516 static const int constr_index[] = {
517 eNR_SHAKE, eNR_SHAKE_RIJ, eNR_SETTLE, eNR_UPDATE, eNR_PCOUPL,
518 eNR_CONSTR_VIR, eNR_CONSTR_V
520 #define NCONSTR_INDEX asize(constr_index)
522 static double pr_av(FILE *log, t_commrec *cr,
523 double fav, double ftot[], const char *title)
525 int i, perc;
526 double dperc, unb;
528 unb = 0;
529 if (fav > 0)
531 fav /= cr->nnodes - cr->npmenodes;
532 fprintf(log, "\n %-26s", title);
533 for (i = 0; (i < cr->nnodes); i++)
535 dperc = (100.0*ftot[i])/fav;
536 unb = std::max(unb, dperc);
537 perc = static_cast<int>(dperc);
538 fprintf(log, "%3d ", perc);
540 if (unb > 0)
542 perc = static_cast<int>(10000.0/unb);
543 fprintf(log, "%6d%%\n\n", perc);
545 else
547 fprintf(log, "\n\n");
550 return unb;
553 void pr_load(FILE *log, t_commrec *cr, t_nrnb nrnb[])
555 int i, j, perc;
556 double dperc, unb, uf, us;
557 double *ftot, fav;
558 double *stot, sav;
559 t_nrnb *av;
561 snew(av, 1);
562 snew(ftot, cr->nnodes);
563 snew(stot, cr->nnodes);
564 init_nrnb(av);
565 for (i = 0; (i < cr->nnodes); i++)
567 add_nrnb(av, av, &(nrnb[i]));
568 /* Cost due to forces */
569 for (j = 0; (j < eNR_NBKERNEL_TOTAL_NR); j++)
571 ftot[i] += nrnb[i].n[j]*cost_nrnb(j);
573 for (j = 0; (j < NFORCE_INDEX); j++)
575 ftot[i] += nrnb[i].n[force_index[j]]*cost_nrnb(force_index[j]);
577 /* Due to shake */
578 for (j = 0; (j < NCONSTR_INDEX); j++)
580 stot[i] += nrnb[i].n[constr_index[j]]*cost_nrnb(constr_index[j]);
583 for (j = 0; (j < eNRNB); j++)
585 av->n[j] = av->n[j]/static_cast<double>(cr->nnodes - cr->npmenodes);
588 fprintf(log, "\nDetailed load balancing info in percentage of average\n");
590 fprintf(log, " Type RANK:");
591 for (i = 0; (i < cr->nnodes); i++)
593 fprintf(log, "%3d ", i);
595 fprintf(log, "Scaling\n");
596 fprintf(log, "---------------------------");
597 for (i = 0; (i < cr->nnodes); i++)
599 fprintf(log, "----");
601 fprintf(log, "-------\n");
603 for (j = 0; (j < eNRNB); j++)
605 unb = 100.0;
606 if (av->n[j] > 0)
608 fprintf(log, " %-26s", nrnb_str(j));
609 for (i = 0; (i < cr->nnodes); i++)
611 dperc = (100.0*nrnb[i].n[j])/av->n[j];
612 unb = std::max(unb, dperc);
613 perc = static_cast<int>(dperc);
614 fprintf(log, "%3d ", perc);
616 if (unb > 0)
618 perc = static_cast<int>(10000.0/unb);
619 fprintf(log, "%6d%%\n", perc);
621 else
623 fprintf(log, "\n");
627 fav = sav = 0;
628 for (i = 0; (i < cr->nnodes); i++)
630 fav += ftot[i];
631 sav += stot[i];
633 uf = pr_av(log, cr, fav, ftot, "Total Force");
634 us = pr_av(log, cr, sav, stot, "Total Constr.");
636 unb = (uf*fav+us*sav)/(fav+sav);
637 if (unb > 0)
639 unb = 10000.0/unb;
640 fprintf(log, "\nTotal Scaling: %.0f%% of max performance\n\n", unb);