Clean up copyrite.*.
[alexxy/gromacs.git] / src / gromacs / gmxlib / nrnb.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <string.h>
40 #include "types/commrec.h"
41 #include "sysstuff.h"
42 #include "gmx_fatal.h"
43 #include "names.h"
44 #include "macros.h"
45 #include "nrnb.h"
46 #include "main.h"
47 #include "smalloc.h"
48
49 typedef struct {
50     const char *name;
51     int         flop;
52 } t_nrnb_data;
53
54
55 static const t_nrnb_data nbdata[eNRNB] = {
56     /* These are re-used for different NB kernels, since there are so many.
57      * The actual number of flops is set dynamically.
58      */
59     { "NB VdW [V&F]",                    1 },
60     { "NB VdW [F]",                      1 },
61     { "NB Elec. [V&F]",                  1 },
62     { "NB Elec. [F]",                    1 },
63     { "NB Elec. [W3,V&F]",               1 },
64     { "NB Elec. [W3,F]",                 1 },
65     { "NB Elec. [W3-W3,V&F]",            1 },
66     { "NB Elec. [W3-W3,F]",              1 },
67     { "NB Elec. [W4,V&F]",               1 },
68     { "NB Elec. [W4,F]",                 1 },
69     { "NB Elec. [W4-W4,V&F]",            1 },
70     { "NB Elec. [W4-W4,F]",              1 },
71     { "NB VdW & Elec. [V&F]",            1 },
72     { "NB VdW & Elec. [F]",              1 },
73     { "NB VdW & Elec. [W3,V&F]",         1 },
74     { "NB VdW & Elec. [W3,F]",           1 },
75     { "NB VdW & Elec. [W3-W3,V&F]",      1 },
76     { "NB VdW & Elec. [W3-W3,F]",        1 },
77     { "NB VdW & Elec. [W4,V&F]",         1 },
78     { "NB VdW & Elec. [W4,F]",           1 },
79     { "NB VdW & Elec. [W4-W4,V&F]",      1 },
80     { "NB VdW & Elec. [W4-W4,F]",        1 },
81
82     { "NB Generic kernel",               1 },
83     { "NB Generic charge grp kernel",    1 },
84     { "NB Generic AdResS kernel",        1 },
85     { "NB Free energy kernel",           1 },
86     { "NB All-vs-all",                   1 },
87     { "NB All-vs-all, GB",               1 },
88
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 CUDA.
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 CUDA, 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.
99      */
100     { "NxN RF Elec. + VdW [F]",         38 }, /* nbnxn kernel LJ+RF, no ener */
101     { "NxN RF Elec. + VdW [V&F]",       54 },
102     { "NxN QSTab Elec. + VdW [F]",      41 }, /* nbnxn kernel LJ+tab, no en */
103     { "NxN QSTab Elec. + VdW [V&F]",    59 },
104     { "NxN Ewald Elec. + VdW [F]",      66 }, /* nbnxn kernel LJ+Ewald, no en */
105     { "NxN Ewald Elec. + VdW [V&F]",   107 },
106     { "NxN VdW [F]",                    33 }, /* nbnxn kernel LJ, no ener */
107     { "NxN VdW [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     { "1,4 nonbonded interactions",     90 },
115     { "Born radii (Still)",             47 },
116     { "Born radii (HCT/OBC)",          183 },
117     { "Born force chain rule",          15 },
118     { "All-vs-All Still radii",          1 },
119     { "All-vs-All HCT/OBC radii",        1 },
120     { "All-vs-All Born chain rule",      1 },
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     { "Mixed Generalized Born stuff",   10 }
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     gmx_bool   bPrint, bPrinted;
203
204     ndays    = dt/(24*3600);
205     dt       = dt-24*3600*ndays;
206     nhours   = dt/3600;
207     dt       = dt-3600*nhours;
208     nmins    = dt/60;
209     dt       = dt-nmins*60;
210     nsecs    = 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 init_nrnb(t_nrnb *nrnb)
255 {
256     int i;
257
258     for (i = 0; (i < eNRNB); i++)
259     {
260         nrnb->n[i] = 0.0;
261     }
262 }
263
264 void cp_nrnb(t_nrnb *dest, t_nrnb *src)
265 {
266     int i;
267
268     for (i = 0; (i < eNRNB); i++)
269     {
270         dest->n[i] = src->n[i];
271     }
272 }
273
274 void add_nrnb(t_nrnb *dest, t_nrnb *s1, t_nrnb *s2)
275 {
276     int i;
277
278     for (i = 0; (i < eNRNB); i++)
279     {
280         dest->n[i] = s1->n[i]+s2->n[i];
281     }
282 }
283
284 void print_nrnb(FILE *out, t_nrnb *nrnb)
285 {
286     int i;
287
288     for (i = 0; (i < eNRNB); i++)
289     {
290         if (nrnb->n[i] > 0)
291         {
292             fprintf(out, " %-26s %10.0f.\n", nbdata[i].name, nrnb->n[i]);
293         }
294     }
295 }
296
297 void _inc_nrnb(t_nrnb *nrnb, int enr, int inc, char *file, int line)
298 {
299     nrnb->n[enr] += inc;
300 #ifdef DEBUG_NRNB
301     printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
302            nbdata[enr].name, enr, inc, file, line);
303 #endif
304 }
305
306 void print_flop(FILE *out, t_nrnb *nrnb, double *nbfs, double *mflop)
307 {
308     int           i;
309     double        mni, frac, tfrac, tflop;
310     const char   *myline = "-----------------------------------------------------------------------------";
311
312     *nbfs = 0.0;
313     for (i = 0; (i < eNR_NBKERNEL_ALLVSALLGB); i++)
314     {
315         if (strstr(nbdata[i].name, "W3-W3") != NULL)
316         {
317             *nbfs += 9e-6*nrnb->n[i];
318         }
319         else if (strstr(nbdata[i].name, "W3") != NULL)
320         {
321             *nbfs += 3e-6*nrnb->n[i];
322         }
323         else if (strstr(nbdata[i].name, "W4-W4") != NULL)
324         {
325             *nbfs += 10e-6*nrnb->n[i];
326         }
327         else if (strstr(nbdata[i].name, "W4") != NULL)
328         {
329             *nbfs += 4e-6*nrnb->n[i];
330         }
331         else
332         {
333             *nbfs += 1e-6*nrnb->n[i];
334         }
335     }
336     tflop = 0;
337     for (i = 0; (i < eNRNB); i++)
338     {
339         tflop += 1e-6*nrnb->n[i]*nbdata[i].flop;
340     }
341
342     if (tflop == 0)
343     {
344         fprintf(out, "No MEGA Flopsen this time\n");
345         return;
346     }
347     if (out)
348     {
349         fprintf(out, "\n\tM E G A - F L O P S   A C C O U N T I N G\n\n");
350     }
351
352     if (out)
353     {
354         fprintf(out, " NB=Group-cutoff nonbonded kernels    NxN=N-by-N cluster Verlet kernels\n");
355         fprintf(out, " RF=Reaction-Field  VdW=Van der Waals  QSTab=quadratic-spline table\n");
356         fprintf(out, " W3=SPC/TIP3p  W4=TIP4p (single or pairs)\n");
357         fprintf(out, " V&F=Potential and force  V=Potential only  F=Force only\n\n");
358
359         fprintf(out, " %-32s %16s %15s  %7s\n",
360                 "Computing:", "M-Number", "M-Flops", "% Flops");
361         fprintf(out, "%s\n", myline);
362     }
363     *mflop = 0.0;
364     tfrac  = 0.0;
365     for (i = 0; (i < eNRNB); i++)
366     {
367         mni     = 1e-6*nrnb->n[i];
368         *mflop += mni*nbdata[i].flop;
369         frac    = 100.0*mni*nbdata[i].flop/tflop;
370         tfrac  += frac;
371         if (out && mni != 0)
372         {
373             fprintf(out, " %-32s %16.6f %15.3f  %6.1f\n",
374                     nbdata[i].name, mni, mni*nbdata[i].flop, frac);
375         }
376     }
377     if (out)
378     {
379         fprintf(out, "%s\n", myline);
380         fprintf(out, " %-32s %16s %15.3f  %6.1f\n",
381                 "Total", "", *mflop, tfrac);
382         fprintf(out, "%s\n\n", myline);
383     }
384 }
385
386 void print_perf(FILE *out, double nodetime, double realtime, int nprocs,
387                 gmx_large_int_t nsteps, real delta_t,
388                 double nbfs, double mflop,
389                 int omp_nth_pp)
390 {
391     real runtime;
392
393     fprintf(out, "\n");
394
395     if (realtime > 0)
396     {
397         fprintf(out, "%12s %12s %12s %10s\n", "", "Core t (s)", "Wall t (s)", "(%)");
398         fprintf(out, "%12s %12.3f %12.3f %10.1f\n", "Time:",
399                 nodetime, realtime, 100.0*nodetime/realtime);
400         /* only print day-hour-sec format if realtime is more than 30 min */
401         if (realtime > 30*60)
402         {
403             fprintf(out, "%12s %12s", "", "");
404             pr_difftime(out, realtime);
405         }
406         if (delta_t > 0)
407         {
408             mflop   = mflop/realtime;
409             runtime = nsteps*delta_t;
410
411             if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
412             {
413                 fprintf(out, "%12s %12s %12s\n",
414                         "", "(ns/day)", "(hour/ns)");
415                 fprintf(out, "%12s %12.3f %12.3f\n", "Performance:",
416                         runtime*24*3.6/realtime, 1000*realtime/(3600*runtime));
417             }
418             else
419             {
420                 fprintf(out, "%12s %12s %12s %12s %12s\n",
421                         "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
422                         "(ns/day)", "(hour/ns)");
423                 fprintf(out, "%12s %12.3f %12.3f %12.3f %12.3f\n", "Performance:",
424                         nbfs/realtime, (mflop > 1000) ? (mflop/1000) : mflop,
425                         runtime*24*3.6/realtime, 1000*realtime/(3600*runtime));
426             }
427         }
428         else
429         {
430             if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
431             {
432                 fprintf(out, "%12s %14s\n",
433                         "", "(steps/hour)");
434                 fprintf(out, "%12s %14.1f\n", "Performance:",
435                         nsteps*3600.0/realtime);
436             }
437             else
438             {
439                 fprintf(out, "%12s %12s %12s %14s\n",
440                         "", "(Mnbf/s)", (mflop > 1000) ? "(GFlops)" : "(MFlops)",
441                         "(steps/hour)");
442                 fprintf(out, "%12s %12.3f %12.3f %14.1f\n", "Performance:",
443                         nbfs/realtime, (mflop > 1000) ? (mflop/1000) : mflop,
444                         nsteps*3600.0/realtime);
445             }
446         }
447     }
448 }
449
450 int cost_nrnb(int enr)
451 {
452     return nbdata[enr].flop;
453 }
454
455 const char *nrnb_str(int enr)
456 {
457     return nbdata[enr].name;
458 }
459
460 static const int    force_index[] = {
461     eNR_BONDS,  eNR_ANGLES,  eNR_PROPER, eNR_IMPROPER,
462     eNR_RB,     eNR_DISRES,  eNR_ORIRES, eNR_POSRES,
463     eNR_FBPOSRES, eNR_NS,
464 };
465 #define NFORCE_INDEX asize(force_index)
466
467 static const int    constr_index[] = {
468     eNR_SHAKE,     eNR_SHAKE_RIJ, eNR_SETTLE,       eNR_UPDATE,       eNR_PCOUPL,
469     eNR_CONSTR_VIR, eNR_CONSTR_V
470 };
471 #define NCONSTR_INDEX asize(constr_index)
472
473 static double pr_av(FILE *log, t_commrec *cr,
474                     double fav, double ftot[], const char *title)
475 {
476     int    i, perc;
477     double dperc, unb;
478
479     unb = 0;
480     if (fav > 0)
481     {
482         fav /= cr->nnodes - cr->npmenodes;
483         fprintf(log, "\n %-26s", title);
484         for (i = 0; (i < cr->nnodes); i++)
485         {
486             dperc = (100.0*ftot[i])/fav;
487             unb   = max(unb, dperc);
488             perc  = dperc;
489             fprintf(log, "%3d ", perc);
490         }
491         if (unb > 0)
492         {
493             perc = 10000.0/unb;
494             fprintf(log, "%6d%%\n\n", perc);
495         }
496         else
497         {
498             fprintf(log, "\n\n");
499         }
500     }
501     return unb;
502 }
503
504 void pr_load(FILE *log, t_commrec *cr, t_nrnb nrnb[])
505 {
506     int     i, j, perc;
507     double  dperc, unb, uf, us;
508     double *ftot, fav;
509     double *stot, sav;
510     t_nrnb *av;
511
512     snew(av, 1);
513     snew(ftot, cr->nnodes);
514     snew(stot, cr->nnodes);
515     init_nrnb(av);
516     for (i = 0; (i < cr->nnodes); i++)
517     {
518         add_nrnb(av, av, &(nrnb[i]));
519         /* Cost due to forces */
520         for (j = 0; (j < eNR_NBKERNEL_ALLVSALLGB); j++)
521         {
522             ftot[i] += nrnb[i].n[j]*cost_nrnb(j);
523         }
524         for (j = 0; (j < NFORCE_INDEX); j++)
525         {
526             ftot[i] += nrnb[i].n[force_index[j]]*cost_nrnb(force_index[j]);
527         }
528         /* Due to shake */
529         for (j = 0; (j < NCONSTR_INDEX); j++)
530         {
531             stot[i] += nrnb[i].n[constr_index[j]]*cost_nrnb(constr_index[j]);
532         }
533     }
534     for (j = 0; (j < eNRNB); j++)
535     {
536         av->n[j] = av->n[j]/(double)(cr->nnodes - cr->npmenodes);
537     }
538
539     fprintf(log, "\nDetailed load balancing info in percentage of average\n");
540
541     fprintf(log, " Type                 NODE:");
542     for (i = 0; (i < cr->nnodes); i++)
543     {
544         fprintf(log, "%3d ", i);
545     }
546     fprintf(log, "Scaling\n");
547     fprintf(log, "---------------------------");
548     for (i = 0; (i < cr->nnodes); i++)
549     {
550         fprintf(log, "----");
551     }
552     fprintf(log, "-------\n");
553
554     for (j = 0; (j < eNRNB); j++)
555     {
556         unb = 100.0;
557         if (av->n[j] > 0)
558         {
559             fprintf(log, " %-26s", nrnb_str(j));
560             for (i = 0; (i < cr->nnodes); i++)
561             {
562                 dperc = (100.0*nrnb[i].n[j])/av->n[j];
563                 unb   = max(unb, dperc);
564                 perc  = dperc;
565                 fprintf(log, "%3d ", perc);
566             }
567             if (unb > 0)
568             {
569                 perc = 10000.0/unb;
570                 fprintf(log, "%6d%%\n", perc);
571             }
572             else
573             {
574                 fprintf(log, "\n");
575             }
576         }
577     }
578     fav = sav = 0;
579     for (i = 0; (i < cr->nnodes); i++)
580     {
581         fav += ftot[i];
582         sav += stot[i];
583     }
584     uf = pr_av(log, cr, fav, ftot, "Total Force");
585     us = pr_av(log, cr, sav, stot, "Total Constr.");
586
587     unb = (uf*fav+us*sav)/(fav+sav);
588     if (unb > 0)
589     {
590         unb = 10000.0/unb;
591         fprintf(log, "\nTotal Scaling: %.0f%% of max performance\n\n", unb);
592     }
593 }