Modernize t_nrnb
[alexxy/gromacs.git] / src / gromacs / gmxlib / nrnb.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
7  * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "nrnb.h"
41
42 #include <cstdlib>
43 #include <cstring>
44
45 typedef struct
46 {
47     const char* name;
48     int         flop;
49 } t_nrnb_data;
50
51
52 static const t_nrnb_data nbdata[eNRNB] = {
53     /* These are re-used for different NB kernels, since there are so many.
54      * The actual number of flops is set dynamically.
55      */
56     { "NB VdW [V&F]", 1 },
57     { "NB VdW [F]", 1 },
58     { "NB Elec. [V&F]", 1 },
59     { "NB Elec. [F]", 1 },
60     { "NB Elec. [W3,V&F]", 1 },
61     { "NB Elec. [W3,F]", 1 },
62     { "NB Elec. [W3-W3,V&F]", 1 },
63     { "NB Elec. [W3-W3,F]", 1 },
64     { "NB Elec. [W4,V&F]", 1 },
65     { "NB Elec. [W4,F]", 1 },
66     { "NB Elec. [W4-W4,V&F]", 1 },
67     { "NB Elec. [W4-W4,F]", 1 },
68     { "NB VdW & Elec. [V&F]", 1 },
69     { "NB VdW & Elec. [F]", 1 },
70     { "NB VdW & Elec. [W3,V&F]", 1 },
71     { "NB VdW & Elec. [W3,F]", 1 },
72     { "NB VdW & Elec. [W3-W3,V&F]", 1 },
73     { "NB VdW & Elec. [W3-W3,F]", 1 },
74     { "NB VdW & Elec. [W4,V&F]", 1 },
75     { "NB VdW & Elec. [W4,F]", 1 },
76     { "NB VdW & Elec. [W4-W4,V&F]", 1 },
77     { "NB VdW & Elec. [W4-W4,F]", 1 },
78
79     { "NB Generic kernel", 1 },
80     { "NB Generic charge grp kernel", 1 },
81     { "NB Free energy kernel", 1 },
82
83     { "Pair Search distance check", 9 }, /* nbnxn pair dist. check */
84     /* nbnxn kernel flops are based on inner-loops without exclusion checks.
85      * Plain Coulomb runs through the RF kernels, except with GPUs.
86      * invsqrt is counted as 6 flops: 1 for _mm_rsqt_ps + 5 for iteration.
87      * The flops are equal for plain-C, x86 SIMD and GPUs, except for:
88      * - plain-C kernel uses one flop more for Coulomb-only (F) than listed
89      * - x86 SIMD LJ geom-comb.rule kernels (fastest) use 2 more flops
90      * - x86 SIMD LJ LB-comb.rule kernels (fast) use 3 (8 for F+E) more flops
91      * - GPU always does exclusions, which requires 2-4 flops, but as invsqrt
92      *   is always counted as 6 flops, this roughly compensates.
93      */
94     { "NxN RF Elec. + LJ [F]", 38 }, /* nbnxn kernel LJ+RF, no ener */
95     { "NxN RF Elec. + LJ [V&F]", 54 },
96     { "NxN QSTab Elec. + LJ [F]", 41 }, /* nbnxn kernel LJ+tab, no en */
97     { "NxN QSTab Elec. + LJ [V&F]", 59 },
98     { "NxN Ewald Elec. + LJ [F]", 66 }, /* nbnxn kernel LJ+Ewald, no en */
99     { "NxN Ewald Elec. + LJ [V&F]", 107 },
100     { "NxN LJ [F]", 33 }, /* nbnxn kernel LJ, no ener */
101     { "NxN LJ [V&F]", 43 },
102     { "NxN RF Electrostatics [F]", 31 }, /* nbnxn kernel RF, no ener */
103     { "NxN RF Electrostatics [V&F]", 36 },
104     { "NxN QSTab Elec. [F]", 34 }, /* nbnxn kernel tab, no ener */
105     { "NxN QSTab Elec. [V&F]", 41 },
106     { "NxN Ewald Elec. [F]", 61 }, /* nbnxn kernel Ewald, no ener */
107     { "NxN Ewald Elec. [V&F]", 84 },
108     /* The switch function flops should be added to the LJ kernels above */
109     { "NxN LJ add F-switch [F]", 12 }, /* extra cost for LJ F-switch */
110     { "NxN LJ add F-switch [V&F]", 22 },
111     { "NxN LJ add P-switch [F]", 27 }, /* extra cost for LJ P-switch */
112     { "NxN LJ add P-switch [V&F]", 20 },
113     { "NxN LJ add LJ Ewald [F]", 36 }, /* extra cost for LJ Ewald */
114     { "NxN LJ add LJ Ewald [V&F]", 33 },
115     { "1,4 nonbonded interactions", 90 },
116     { "Calc Weights", 36 },
117     { "Spread Q", 6 },
118     { "Spread Q Bspline", 2 },
119     { "Gather F", 23 },
120     { "Gather F Bspline", 6 },
121     { "3D-FFT", 8 },
122     { "Convolution", 4 },
123     { "Solve PME", 64 },
124     { "NS-Pairs", 21 },
125     { "Reset In Box", 3 },
126     { "Shift-X", 6 },
127     { "CG-CoM", 3 },
128     { "Sum Forces", 1 },
129     { "Bonds", 59 },
130     { "G96Bonds", 44 },
131     { "FENE Bonds", 58 },
132     { "Tab. Bonds", 62 },
133     { "Restraint Potential", 86 },
134     { "Linear Angles", 57 },
135     { "Angles", 168 },
136     { "G96Angles", 150 },
137     { "Quartic Angles", 160 },
138     { "Tab. Angles", 169 },
139     { "Propers", 229 },
140     { "Impropers", 208 },
141     { "RB-Dihedrals", 247 },
142     { "Four. Dihedrals", 247 },
143     { "Tab. Dihedrals", 227 },
144     { "Dist. Restr.", 200 },
145     { "Orient. Restr.", 200 },
146     { "Dihedral Restr.", 200 },
147     { "Pos. Restr.", 50 },
148     { "Flat-bottom posres", 50 },
149     { "Angle Restr.", 191 },
150     { "Angle Restr. Z", 164 },
151     { "Morse Potent.", 83 },
152     { "Cubic Bonds", 54 },
153     { "Walls", 31 },
154     { "Polarization", 59 },
155     { "Anharmonic Polarization", 72 },
156     { "Water Pol.", 62 },
157     { "Thole Pol.", 296 },
158     { "Virial", 18 },
159     { "Update", 31 },
160     { "Ext.ens. Update", 54 },
161     { "Stop-CM", 10 },
162     { "P-Coupling", 6 },
163     { "Calc-Ekin", 27 },
164     { "Lincs", 60 },
165     { "Lincs-Mat", 4 },
166     { "Shake", 30 },
167     { "Constraint-V", 9 },
168     { "Shake-Init", 10 },
169     { "Constraint-Vir", 24 },
170     { "Settle", 370 },
171     { "Virtual Site 1", 1 },
172     { "Virtual Site 2", 23 },
173     { "Virtual Site 2fd", 63 },
174     { "Virtual Site 3", 37 },
175     { "Virtual Site 3fd", 95 },
176     { "Virtual Site 3fad", 176 },
177     { "Virtual Site 3out", 87 },
178     { "Virtual Site 4fd", 110 },
179     { "Virtual Site 4fdn", 254 },
180     { "Virtual Site N", 15 },
181     { "CMAP", 1700 }, // Estimate!
182     { "Urey-Bradley", 183 },
183     { "Cross-Bond-Bond", 163 },
184     { "Cross-Bond-Angle", 163 }
185 };
186
187 static void pr_two(FILE* out, int c, int i)
188 {
189     if (i < 10)
190     {
191         fprintf(out, "%c0%1d", c, i);
192     }
193     else
194     {
195         fprintf(out, "%c%2d", c, i);
196     }
197 }
198
199 static void pr_difftime(FILE* out, double dt)
200 {
201     int  ndays, nhours, nmins, nsecs;
202     bool bPrint, bPrinted;
203
204     ndays    = static_cast<int>(dt / (24 * 3600));
205     dt       = dt - 24 * 3600 * ndays;
206     nhours   = static_cast<int>(dt / 3600);
207     dt       = dt - 3600 * nhours;
208     nmins    = static_cast<int>(dt / 60);
209     dt       = dt - nmins * 60;
210     nsecs    = static_cast<int>(dt);
211     bPrint   = (ndays > 0);
212     bPrinted = bPrint;
213     if (bPrint)
214     {
215         fprintf(out, "%d", ndays);
216     }
217     bPrint = bPrint || (nhours > 0);
218     if (bPrint)
219     {
220         if (bPrinted)
221         {
222             pr_two(out, 'd', nhours);
223         }
224         else
225         {
226             fprintf(out, "%d", nhours);
227         }
228     }
229     bPrinted = bPrinted || bPrint;
230     bPrint   = bPrint || (nmins > 0);
231     if (bPrint)
232     {
233         if (bPrinted)
234         {
235             pr_two(out, 'h', nmins);
236         }
237         else
238         {
239             fprintf(out, "%d", nmins);
240         }
241     }
242     bPrinted = bPrinted || bPrint;
243     if (bPrinted)
244     {
245         pr_two(out, ':', nsecs);
246     }
247     else
248     {
249         fprintf(out, "%ds", nsecs);
250     }
251     fprintf(out, "\n");
252 }
253
254 void clear_nrnb(t_nrnb* nrnb)
255 {
256     for (int i = 0; (i < eNRNB); i++)
257     {
258         nrnb->n[i] = 0.0;
259     }
260 }
261
262 void print_nrnb(FILE* out, t_nrnb* nrnb)
263 {
264     for (int i = 0; (i < eNRNB); i++)
265     {
266         if (nrnb->n[i] > 0)
267         {
268             fprintf(out, " %-26s %10.0f.\n", nbdata[i].name, nrnb->n[i]);
269         }
270     }
271 }
272
273 /* Returns in enr is the index of a full nbnxn VdW kernel */
274 static bool nrnb_is_nbnxn_vdw_kernel(int enr)
275 {
276     return (enr >= eNR_NBNXN_LJ_RF && enr <= eNR_NBNXN_LJ_E);
277 }
278
279 /* Returns in enr is the index of an nbnxn kernel addition (LJ modification) */
280 static bool nrnb_is_nbnxn_kernel_addition(int enr)
281 {
282     return (enr >= eNR_NBNXN_ADD_LJ_FSW && enr <= eNR_NBNXN_ADD_LJ_EWALD_E);
283 }
284
285 void atomicNrnbIncrement(t_nrnb* nrnb, int index, int increment)
286 {
287 #pragma omp atomic
288     nrnb->n[index] += increment;
289 }
290
291 void print_flop(FILE* out, t_nrnb* nrnb, double* nbfs, double* mflop)
292 {
293     double      mni, frac, tfrac, tflop;
294     const char* myline =
295             "-----------------------------------------------------------------------------";
296
297     *nbfs = 0.0;
298     for (int i = 0; (i < eNR_NBKERNEL_TOTAL_NR); i++)
299     {
300         if (std::strstr(nbdata[i].name, "W3-W3") != nullptr)
301         {
302             *nbfs += 9e-6 * nrnb->n[i];
303         }
304         else if (std::strstr(nbdata[i].name, "W3") != nullptr)
305         {
306             *nbfs += 3e-6 * nrnb->n[i];
307         }
308         else if (std::strstr(nbdata[i].name, "W4-W4") != nullptr)
309         {
310             *nbfs += 10e-6 * nrnb->n[i];
311         }
312         else if (std::strstr(nbdata[i].name, "W4") != nullptr)
313         {
314             *nbfs += 4e-6 * nrnb->n[i];
315         }
316         else
317         {
318             *nbfs += 1e-6 * nrnb->n[i];
319         }
320     }
321     tflop = 0;
322     for (int i = 0; (i < eNRNB); i++)
323     {
324         tflop += 1e-6 * nrnb->n[i] * nbdata[i].flop;
325     }
326
327     if (tflop == 0)
328     {
329         fprintf(out, "No MEGA Flopsen this time\n");
330         return;
331     }
332     if (out)
333     {
334         fprintf(out, "\n\tM E G A - F L O P S   A C C O U N T I N G\n\n");
335     }
336
337     if (out)
338     {
339         fprintf(out, " NB=Group-cutoff nonbonded kernels    NxN=N-by-N cluster Verlet kernels\n");
340         fprintf(out, " RF=Reaction-Field  VdW=Van der Waals  QSTab=quadratic-spline table\n");
341         fprintf(out, " W3=SPC/TIP3p  W4=TIP4p (single or pairs)\n");
342         fprintf(out, " V&F=Potential and force  V=Potential only  F=Force only\n\n");
343
344         fprintf(out, " %-32s %16s %15s  %7s\n", "Computing:", "M-Number", "M-Flops", "% Flops");
345         fprintf(out, "%s\n", myline);
346     }
347     *mflop = 0.0;
348     tfrac  = 0.0;
349     for (int i = 0; (i < eNRNB); i++)
350     {
351         mni = 1e-6 * nrnb->n[i];
352         /* Skip empty entries and nbnxn additional flops,
353          * which have been added to the kernel entry.
354          */
355         if (mni > 0 && !nrnb_is_nbnxn_kernel_addition(i))
356         {
357             int flop;
358
359             flop = nbdata[i].flop;
360             if (nrnb_is_nbnxn_vdw_kernel(i))
361             {
362                 /* Possibly add the cost of an LJ switch/Ewald function */
363                 for (int j = eNR_NBNXN_ADD_LJ_FSW; j <= eNR_NBNXN_ADD_LJ_EWALD; j += 2)
364                 {
365                     int e_kernel_add;
366
367                     /* Select the force or energy flop count */
368                     e_kernel_add = j + ((i - eNR_NBNXN_LJ_RF) % 2);
369
370                     if (nrnb->n[e_kernel_add] > 0)
371                     {
372                         flop += nbdata[e_kernel_add].flop;
373                     }
374                 }
375             }
376             *mflop += mni * flop;
377             frac = 100.0 * mni * flop / tflop;
378             tfrac += frac;
379             if (out != nullptr)
380             {
381                 fprintf(out, " %-32s %16.6f %15.3f  %6.1f\n", nbdata[i].name, mni, mni * flop, frac);
382             }
383         }
384     }
385     if (out)
386     {
387         fprintf(out, "%s\n", myline);
388         fprintf(out, " %-32s %16s %15.3f  %6.1f\n", "Total", "", *mflop, tfrac);
389         fprintf(out, "%s\n\n", myline);
390
391         if (nrnb->n[eNR_NBKERNEL_GENERIC] > 0)
392         {
393             fprintf(out,
394                     "WARNING: Using the slow generic C kernel. This is fine if you are\n"
395                     "comparing different implementations or MD software. Routine\n"
396                     "simulations should use a different non-bonded setup for much better\n"
397                     "performance.\n\n");
398         }
399     }
400 }
401
402 void print_perf(FILE*   out,
403                 double  time_per_thread,
404                 double  time_per_node,
405                 int64_t nsteps,
406                 double  delta_t,
407                 double  nbfs,
408                 double  mflop)
409 {
410     double wallclocktime;
411
412     fprintf(out, "\n");
413
414     if (time_per_node > 0)
415     {
416         fprintf(out, "%12s %12s %12s %10s\n", "", "Core t (s)", "Wall t (s)", "(%)");
417         fprintf(out, "%12s %12.3f %12.3f %10.1f\n", "Time:", time_per_thread, time_per_node, 100.0 * time_per_thread / time_per_node);
418         /* only print day-hour-sec format if time_per_node is more than 30 min */
419         if (time_per_node > 30 * 60)
420         {
421             fprintf(out, "%12s %12s", "", "");
422             pr_difftime(out, time_per_node);
423         }
424         if (delta_t > 0)
425         {
426             mflop         = mflop / time_per_node;
427             wallclocktime = nsteps * delta_t;
428
429             if (getenv("GMX_DETAILED_PERF_STATS") == nullptr)
430             {
431                 fprintf(out, "%12s %12s %12s\n", "", "(ns/day)", "(hour/ns)");
432                 fprintf(out,
433                         "%12s %12.3f %12.3f\n",
434                         "Performance:",
435                         wallclocktime * 24 * 3.6 / time_per_node,
436                         1000 * time_per_node / (3600 * wallclocktime));
437             }
438             else
439             {
440                 fprintf(out,
441                         "%12s %12s %12s %12s %12s\n",
442                         "",
443                         "(Mnbf/s)",
444                         (mflop > 1000) ? "(GFlops)" : "(MFlops)",
445                         "(ns/day)",
446                         "(hour/ns)");
447                 fprintf(out,
448                         "%12s %12.3f %12.3f %12.3f %12.3f\n",
449                         "Performance:",
450                         nbfs / time_per_node,
451                         (mflop > 1000) ? (mflop / 1000) : mflop,
452                         wallclocktime * 24 * 3.6 / time_per_node,
453                         1000 * time_per_node / (3600 * wallclocktime));
454             }
455         }
456         else
457         {
458             if (getenv("GMX_DETAILED_PERF_STATS") == nullptr)
459             {
460                 fprintf(out, "%12s %14s\n", "", "(steps/hour)");
461                 fprintf(out, "%12s %14.1f\n", "Performance:", nsteps * 3600.0 / time_per_node);
462             }
463             else
464             {
465                 fprintf(out,
466                         "%12s %12s %12s %14s\n",
467                         "",
468                         "(Mnbf/s)",
469                         (mflop > 1000) ? "(GFlops)" : "(MFlops)",
470                         "(steps/hour)");
471                 fprintf(out,
472                         "%12s %12.3f %12.3f %14.1f\n",
473                         "Performance:",
474                         nbfs / time_per_node,
475                         (mflop > 1000) ? (mflop / 1000) : mflop,
476                         nsteps * 3600.0 / time_per_node);
477             }
478         }
479     }
480 }
481
482 int cost_nrnb(int enr)
483 {
484     return nbdata[enr].flop;
485 }
486
487 const char* nrnb_str(int enr)
488 {
489     return nbdata[enr].name;
490 }