Fix remaining copyright headers
[alexxy/gromacs.git] / src / gromacs / gmxlib / nrnb.c
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, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <string.h>
42 #include "types/commrec.h"
43 #include "sysstuff.h"
44 #include "gmx_fatal.h"
45 #include "names.h"
46 #include "macros.h"
47 #include "nrnb.h"
48 #include "main.h"
49 #include "smalloc.h"
50
51 typedef struct {
52     const char *name;
53     int         flop;
54 } t_nrnb_data;
55
56
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.
60      */
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 },
83
84     { "NB Generic kernel",               1 },
85     { "NB Generic charge grp kernel",    1 },
86     { "NB Generic AdResS kernel",        1 },
87     { "NB Free energy kernel",           1 },
88     { "NB All-vs-all",                   1 },
89     { "NB All-vs-all, GB",               1 },
90
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 CUDA.
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 CUDA, 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.
101      */
102     { "NxN RF Elec. + VdW [F]",         38 }, /* nbnxn kernel LJ+RF, no ener */
103     { "NxN RF Elec. + VdW [V&F]",       54 },
104     { "NxN QSTab Elec. + VdW [F]",      41 }, /* nbnxn kernel LJ+tab, no en */
105     { "NxN QSTab Elec. + VdW [V&F]",    59 },
106     { "NxN Ewald Elec. + VdW [F]",      66 }, /* nbnxn kernel LJ+Ewald, no en */
107     { "NxN Ewald Elec. + VdW [V&F]",   107 },
108     { "NxN VdW [F]",                    33 }, /* nbnxn kernel LJ, no ener */
109     { "NxN VdW [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     { "1,4 nonbonded interactions",     90 },
117     { "Born radii (Still)",             47 },
118     { "Born radii (HCT/OBC)",          183 },
119     { "Born force chain rule",          15 },
120     { "All-vs-All Still radii",          1 },
121     { "All-vs-All HCT/OBC radii",        1 },
122     { "All-vs-All Born chain rule",      1 },
123     { "Calc Weights",                   36 },
124     { "Spread Q",                        6 },
125     { "Spread Q Bspline",                2 },
126     { "Gather F",                      23  },
127     { "Gather F Bspline",              6   },
128     { "3D-FFT",                        8   },
129     { "Convolution",                   4   },
130     { "Solve PME",                     64  },
131     { "NS-Pairs",                      21  },
132     { "Reset In Box",                  3   },
133     { "Shift-X",                       6   },
134     { "CG-CoM",                        3   },
135     { "Sum Forces",                    1   },
136     { "Bonds",                         59  },
137     { "G96Bonds",                      44  },
138     { "FENE Bonds",                    58  },
139     { "Tab. Bonds",                    62  },
140     { "Restraint Potential",           86  },
141     { "Linear Angles",                 57  },
142     { "Angles",                        168 },
143     { "G96Angles",                     150 },
144     { "Quartic Angles",                160 },
145     { "Tab. Angles",                   169 },
146     { "Propers",                       229 },
147     { "Impropers",                     208 },
148     { "RB-Dihedrals",                  247 },
149     { "Four. Dihedrals",               247 },
150     { "Tab. Dihedrals",                227 },
151     { "Dist. Restr.",                  200 },
152     { "Orient. Restr.",                200 },
153     { "Dihedral Restr.",               200 },
154     { "Pos. Restr.",                   50  },
155     { "Flat-bottom posres",            50  },
156     { "Angle Restr.",                  191 },
157     { "Angle Restr. Z",                164 },
158     { "Morse Potent.",                 83  },
159     { "Cubic Bonds",                   54  },
160     { "Walls",                         31  },
161     { "Polarization",                  59  },
162     { "Anharmonic Polarization",       72  },
163     { "Water Pol.",                    62  },
164     { "Thole Pol.",                    296 },
165     { "Virial",                        18  },
166     { "Update",                        31  },
167     { "Ext.ens. Update",               54  },
168     { "Stop-CM",                       10  },
169     { "P-Coupling",                    6   },
170     { "Calc-Ekin",                     27  },
171     { "Lincs",                         60  },
172     { "Lincs-Mat",                     4   },
173     { "Shake",                         30  },
174     { "Constraint-V",                   8  },
175     { "Shake-Init",                    10  },
176     { "Constraint-Vir",                24  },
177     { "Settle",                        323 },
178     { "Virtual Site 2",                23  },
179     { "Virtual Site 3",                37  },
180     { "Virtual Site 3fd",              95  },
181     { "Virtual Site 3fad",             176 },
182     { "Virtual Site 3out",             87  },
183     { "Virtual Site 4fd",              110 },
184     { "Virtual Site 4fdn",             254 },
185     { "Virtual Site N",                 15 },
186     { "Mixed Generalized Born stuff",   10 }
187 };
188
189 static void pr_two(FILE *out, int c, int i)
190 {
191     if (i < 10)
192     {
193         fprintf(out, "%c0%1d", c, i);
194     }
195     else
196     {
197         fprintf(out, "%c%2d", c, i);
198     }
199 }
200
201 static void pr_difftime(FILE *out, double dt)
202 {
203     int        ndays, nhours, nmins, nsecs;
204     gmx_bool   bPrint, bPrinted;
205
206     ndays    = dt/(24*3600);
207     dt       = dt-24*3600*ndays;
208     nhours   = dt/3600;
209     dt       = dt-3600*nhours;
210     nmins    = dt/60;
211     dt       = dt-nmins*60;
212     nsecs    = dt;
213     bPrint   = (ndays > 0);
214     bPrinted = bPrint;
215     if (bPrint)
216     {
217         fprintf(out, "%d", ndays);
218     }
219     bPrint = bPrint || (nhours > 0);
220     if (bPrint)
221     {
222         if (bPrinted)
223         {
224             pr_two(out, 'd', nhours);
225         }
226         else
227         {
228             fprintf(out, "%d", nhours);
229         }
230     }
231     bPrinted = bPrinted || bPrint;
232     bPrint   = bPrint || (nmins > 0);
233     if (bPrint)
234     {
235         if (bPrinted)
236         {
237             pr_two(out, 'h', nmins);
238         }
239         else
240         {
241             fprintf(out, "%d", nmins);
242         }
243     }
244     bPrinted = bPrinted || bPrint;
245     if (bPrinted)
246     {
247         pr_two(out, ':', nsecs);
248     }
249     else
250     {
251         fprintf(out, "%ds", nsecs);
252     }
253     fprintf(out, "\n");
254 }
255
256 void init_nrnb(t_nrnb *nrnb)
257 {
258     int i;
259
260     for (i = 0; (i < eNRNB); i++)
261     {
262         nrnb->n[i] = 0.0;
263     }
264 }
265
266 void cp_nrnb(t_nrnb *dest, t_nrnb *src)
267 {
268     int i;
269
270     for (i = 0; (i < eNRNB); i++)
271     {
272         dest->n[i] = src->n[i];
273     }
274 }
275
276 void add_nrnb(t_nrnb *dest, t_nrnb *s1, t_nrnb *s2)
277 {
278     int i;
279
280     for (i = 0; (i < eNRNB); i++)
281     {
282         dest->n[i] = s1->n[i]+s2->n[i];
283     }
284 }
285
286 void print_nrnb(FILE *out, t_nrnb *nrnb)
287 {
288     int i;
289
290     for (i = 0; (i < eNRNB); i++)
291     {
292         if (nrnb->n[i] > 0)
293         {
294             fprintf(out, " %-26s %10.0f.\n", nbdata[i].name, nrnb->n[i]);
295         }
296     }
297 }
298
299 void _inc_nrnb(t_nrnb *nrnb, int enr, int inc, char gmx_unused *file, int gmx_unused line)
300 {
301     nrnb->n[enr] += inc;
302 #ifdef DEBUG_NRNB
303     printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
304            nbdata[enr].name, enr, inc, file, line);
305 #endif
306 }
307
308 void print_flop(FILE *out, t_nrnb *nrnb, double *nbfs, double *mflop)
309 {
310     int           i;
311     double        mni, frac, tfrac, tflop;
312     const char   *myline = "-----------------------------------------------------------------------------";
313
314     *nbfs = 0.0;
315     for (i = 0; (i < eNR_NBKERNEL_ALLVSALLGB); i++)
316     {
317         if (strstr(nbdata[i].name, "W3-W3") != NULL)
318         {
319             *nbfs += 9e-6*nrnb->n[i];
320         }
321         else if (strstr(nbdata[i].name, "W3") != NULL)
322         {
323             *nbfs += 3e-6*nrnb->n[i];
324         }
325         else if (strstr(nbdata[i].name, "W4-W4") != NULL)
326         {
327             *nbfs += 10e-6*nrnb->n[i];
328         }
329         else if (strstr(nbdata[i].name, "W4") != NULL)
330         {
331             *nbfs += 4e-6*nrnb->n[i];
332         }
333         else
334         {
335             *nbfs += 1e-6*nrnb->n[i];
336         }
337     }
338     tflop = 0;
339     for (i = 0; (i < eNRNB); i++)
340     {
341         tflop += 1e-6*nrnb->n[i]*nbdata[i].flop;
342     }
343
344     if (tflop == 0)
345     {
346         fprintf(out, "No MEGA Flopsen this time\n");
347         return;
348     }
349     if (out)
350     {
351         fprintf(out, "\n\tM E G A - F L O P S   A C C O U N T I N G\n\n");
352     }
353
354     if (out)
355     {
356         fprintf(out, " NB=Group-cutoff nonbonded kernels    NxN=N-by-N cluster Verlet kernels\n");
357         fprintf(out, " RF=Reaction-Field  VdW=Van der Waals  QSTab=quadratic-spline table\n");
358         fprintf(out, " W3=SPC/TIP3p  W4=TIP4p (single or pairs)\n");
359         fprintf(out, " V&F=Potential and force  V=Potential only  F=Force only\n\n");
360
361         fprintf(out, " %-32s %16s %15s  %7s\n",
362                 "Computing:", "M-Number", "M-Flops", "% Flops");
363         fprintf(out, "%s\n", myline);
364     }
365     *mflop = 0.0;
366     tfrac  = 0.0;
367     for (i = 0; (i < eNRNB); i++)
368     {
369         mni     = 1e-6*nrnb->n[i];
370         *mflop += mni*nbdata[i].flop;
371         frac    = 100.0*mni*nbdata[i].flop/tflop;
372         tfrac  += frac;
373         if (out && mni != 0)
374         {
375             fprintf(out, " %-32s %16.6f %15.3f  %6.1f\n",
376                     nbdata[i].name, mni, mni*nbdata[i].flop, frac);
377         }
378     }
379     if (out)
380     {
381         fprintf(out, "%s\n", myline);
382         fprintf(out, " %-32s %16s %15.3f  %6.1f\n",
383                 "Total", "", *mflop, tfrac);
384         fprintf(out, "%s\n\n", myline);
385     }
386 }
387
388 void print_perf(FILE *out, double time_per_thread, double time_per_node,
389                 gmx_int64_t nsteps, real delta_t,
390                 double nbfs, double mflop)
391 {
392     real wallclocktime;
393
394     fprintf(out, "\n");
395
396     if (time_per_node > 0)
397     {
398         fprintf(out, "%12s %12s %12s %10s\n", "", "Core t (s)", "Wall t (s)", "(%)");
399         fprintf(out, "%12s %12.3f %12.3f %10.1f\n", "Time:",
400                 time_per_thread, time_per_node, 100.0*time_per_thread/time_per_node);
401         /* only print day-hour-sec format if time_per_node is more than 30 min */
402         if (time_per_node > 30*60)
403         {
404             fprintf(out, "%12s %12s", "", "");
405             pr_difftime(out, time_per_node);
406         }
407         if (delta_t > 0)
408         {
409             mflop         = mflop/time_per_node;
410             wallclocktime = nsteps*delta_t;
411
412             if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
413             {
414                 fprintf(out, "%12s %12s %12s\n",
415                         "", "(ns/day)", "(hour/ns)");
416                 fprintf(out, "%12s %12.3f %12.3f\n", "Performance:",
417                         wallclocktime*24*3.6/time_per_node, 1000*time_per_node/(3600*wallclocktime));
418             }
419             else
420             {
421                 fprintf(out, "%12s %12s %12s %12s %12s\n",
422                         "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
423                         "(ns/day)", "(hour/ns)");
424                 fprintf(out, "%12s %12.3f %12.3f %12.3f %12.3f\n", "Performance:",
425                         nbfs/time_per_node, (mflop > 1000) ? (mflop/1000) : mflop,
426                         wallclocktime*24*3.6/time_per_node, 1000*time_per_node/(3600*wallclocktime));
427             }
428         }
429         else
430         {
431             if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
432             {
433                 fprintf(out, "%12s %14s\n",
434                         "", "(steps/hour)");
435                 fprintf(out, "%12s %14.1f\n", "Performance:",
436                         nsteps*3600.0/time_per_node);
437             }
438             else
439             {
440                 fprintf(out, "%12s %12s %12s %14s\n",
441                         "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
442                         "(steps/hour)");
443                 fprintf(out, "%12s %12.3f %12.3f %14.1f\n", "Performance:",
444                         nbfs/time_per_node, (mflop > 1000) ? (mflop/1000) : mflop,
445                         nsteps*3600.0/time_per_node);
446             }
447         }
448     }
449 }
450
451 int cost_nrnb(int enr)
452 {
453     return nbdata[enr].flop;
454 }
455
456 const char *nrnb_str(int enr)
457 {
458     return nbdata[enr].name;
459 }
460
461 static const int    force_index[] = {
462     eNR_BONDS,  eNR_ANGLES,  eNR_PROPER, eNR_IMPROPER,
463     eNR_RB,     eNR_DISRES,  eNR_ORIRES, eNR_POSRES,
464     eNR_FBPOSRES, eNR_NS,
465 };
466 #define NFORCE_INDEX asize(force_index)
467
468 static const int    constr_index[] = {
469     eNR_SHAKE,     eNR_SHAKE_RIJ, eNR_SETTLE,       eNR_UPDATE,       eNR_PCOUPL,
470     eNR_CONSTR_VIR, eNR_CONSTR_V
471 };
472 #define NCONSTR_INDEX asize(constr_index)
473
474 static double pr_av(FILE *log, t_commrec *cr,
475                     double fav, double ftot[], const char *title)
476 {
477     int    i, perc;
478     double dperc, unb;
479
480     unb = 0;
481     if (fav > 0)
482     {
483         fav /= cr->nnodes - cr->npmenodes;
484         fprintf(log, "\n %-26s", title);
485         for (i = 0; (i < cr->nnodes); i++)
486         {
487             dperc = (100.0*ftot[i])/fav;
488             unb   = max(unb, dperc);
489             perc  = dperc;
490             fprintf(log, "%3d ", perc);
491         }
492         if (unb > 0)
493         {
494             perc = 10000.0/unb;
495             fprintf(log, "%6d%%\n\n", perc);
496         }
497         else
498         {
499             fprintf(log, "\n\n");
500         }
501     }
502     return unb;
503 }
504
505 void pr_load(FILE *log, t_commrec *cr, t_nrnb nrnb[])
506 {
507     int     i, j, perc;
508     double  dperc, unb, uf, us;
509     double *ftot, fav;
510     double *stot, sav;
511     t_nrnb *av;
512
513     snew(av, 1);
514     snew(ftot, cr->nnodes);
515     snew(stot, cr->nnodes);
516     init_nrnb(av);
517     for (i = 0; (i < cr->nnodes); i++)
518     {
519         add_nrnb(av, av, &(nrnb[i]));
520         /* Cost due to forces */
521         for (j = 0; (j < eNR_NBKERNEL_ALLVSALLGB); j++)
522         {
523             ftot[i] += nrnb[i].n[j]*cost_nrnb(j);
524         }
525         for (j = 0; (j < NFORCE_INDEX); j++)
526         {
527             ftot[i] += nrnb[i].n[force_index[j]]*cost_nrnb(force_index[j]);
528         }
529         /* Due to shake */
530         for (j = 0; (j < NCONSTR_INDEX); j++)
531         {
532             stot[i] += nrnb[i].n[constr_index[j]]*cost_nrnb(constr_index[j]);
533         }
534     }
535     for (j = 0; (j < eNRNB); j++)
536     {
537         av->n[j] = av->n[j]/(double)(cr->nnodes - cr->npmenodes);
538     }
539
540     fprintf(log, "\nDetailed load balancing info in percentage of average\n");
541
542     fprintf(log, " Type                 NODE:");
543     for (i = 0; (i < cr->nnodes); i++)
544     {
545         fprintf(log, "%3d ", i);
546     }
547     fprintf(log, "Scaling\n");
548     fprintf(log, "---------------------------");
549     for (i = 0; (i < cr->nnodes); i++)
550     {
551         fprintf(log, "----");
552     }
553     fprintf(log, "-------\n");
554
555     for (j = 0; (j < eNRNB); j++)
556     {
557         unb = 100.0;
558         if (av->n[j] > 0)
559         {
560             fprintf(log, " %-26s", nrnb_str(j));
561             for (i = 0; (i < cr->nnodes); i++)
562             {
563                 dperc = (100.0*nrnb[i].n[j])/av->n[j];
564                 unb   = max(unb, dperc);
565                 perc  = dperc;
566                 fprintf(log, "%3d ", perc);
567             }
568             if (unb > 0)
569             {
570                 perc = 10000.0/unb;
571                 fprintf(log, "%6d%%\n", perc);
572             }
573             else
574             {
575                 fprintf(log, "\n");
576             }
577         }
578     }
579     fav = sav = 0;
580     for (i = 0; (i < cr->nnodes); i++)
581     {
582         fav += ftot[i];
583         sav += stot[i];
584     }
585     uf = pr_av(log, cr, fav, ftot, "Total Force");
586     us = pr_av(log, cr, sav, stot, "Total Constr.");
587
588     unb = (uf*fav+us*sav)/(fav+sav);
589     if (unb > 0)
590     {
591         unb = 10000.0/unb;
592         fprintf(log, "\nTotal Scaling: %.0f%% of max performance\n\n", unb);
593     }
594 }