09c56faea7436048530b4c7a4a318c475679bf8b
[alexxy/gromacs.git] / src / 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 #include "copyrite.h"
49
50
51
52
53
54 typedef struct {
55   const char *name;
56   int  flop;
57 } t_nrnb_data;
58
59
60 static const t_nrnb_data nbdata[eNRNB] = {
61     /* These are re-used for different NB kernels, since there are so many.
62      * The actual number of flops is set dynamically.
63      */
64     { "NB VdW [V&F]",                    1 },
65     { "NB VdW [F]",                      1 },
66     { "NB Elec. [V&F]",                  1 },
67     { "NB Elec. [F]",                    1 },
68     { "NB Elec. [W3,V&F]",               1 },
69     { "NB Elec. [W3,F]",                 1 },
70     { "NB Elec. [W3-W3,V&F]",            1 },
71     { "NB Elec. [W3-W3,F]",              1 },
72     { "NB Elec. [W4,V&F]",               1 },
73     { "NB Elec. [W4,F]",                 1 },
74     { "NB Elec. [W4-W4,V&F]",            1 },
75     { "NB Elec. [W4-W4,F]",              1 },
76     { "NB VdW & Elec. [V&F]",            1 },
77     { "NB VdW & Elec. [F]",              1 },
78     { "NB VdW & Elec. [W3,V&F]",         1 },
79     { "NB VdW & Elec. [W3,F]",           1 },
80     { "NB VdW & Elec. [W3-W3,V&F]",      1 },
81     { "NB VdW & Elec. [W3-W3,F]",        1 },
82     { "NB VdW & Elec. [W4,V&F]",         1 },
83     { "NB VdW & Elec. [W4,F]",           1 },
84     { "NB VdW & Elec. [W4-W4,V&F]",      1 },
85     { "NB VdW & Elec. [W4-W4,F]",        1 },
86     
87     { "NB Generic kernel",               1 },
88     { "NB Free energy kernel",           1 },
89     { "NB All-vs-all",                   1 },
90     { "NB All-vs-all, GB",               1 },
91
92     { "Pair Search distance check",      9 }, /* nbnxn pair dist. check */
93     /* nbnxn kernel flops are based on inner-loops without exclusion checks.
94      * Plain Coulomb runs through the RF kernels, except with CUDA.
95      * invsqrt is counted as 6 flops: 1 for _mm_rsqt_ps + 5 for iteration.
96      * The flops are equal for plain-C, x86 SIMD and CUDA, except for:
97      * - plain-C kernel uses one flop more for Coulomb-only (F) than listed
98      * - x86 SIMD LJ geom-comb.rule kernels (fastest) use 2 more flops
99      * - x86 SIMD LJ LB-comb.rule kernels (fast) use 3 (8 for F+E) more flops
100      * - GPU always does exclusions, which requires 2-4 flops, but as invsqrt
101      *   is always counted as 6 flops, this roughly compensates.
102      */
103     { "NxN RF Elec. + VdW [F]",         38 }, /* nbnxn kernel LJ+RF, no ener */
104     { "NxN RF Elec. + VdW [V&F]",       54 },
105     { "NxN QSTab Elec. + VdW [F]",      41 }, /* nbnxn kernel LJ+tab, no en */
106     { "NxN QSTab Elec. + VdW [V&F]",    59 },
107     { "NxN Ewald Elec. + VdW [F]",      66 }, /* nbnxn kernel LJ+Ewald, no en */
108     { "NxN Ewald Elec. + VdW [V&F]",   107 },
109     { "NxN VdW [F]",                    33 }, /* nbnxn kernel LJ, no ener */
110     { "NxN VdW [V&F]",                  43 },
111     { "NxN RF Electrostatics [F]",      31 }, /* nbnxn kernel RF, no ener */
112     { "NxN RF Electrostatics [V&F]",    36 },
113     { "NxN QSTab Elec. [F]",            34 }, /* nbnxn kernel tab, no ener */
114     { "NxN QSTab Elec. [V&F]",          41 },
115     { "NxN Ewald Elec. [F]",            61 }, /* nbnxn kernel Ewald, no ener */
116     { "NxN Ewald Elec. [V&F]",          84 },
117     { "1,4 nonbonded interactions",     90 },
118     { "Born radii (Still)",             47 },
119     { "Born radii (HCT/OBC)",          183 },
120     { "Born force chain rule",          15 },
121     { "All-vs-All Still radii",          1 },
122     { "All-vs-All HCT/OBC radii",        1 },
123     { "All-vs-All Born chain rule",      1 },
124     { "Calc Weights",                   36 },
125     { "Spread Q",                        6 },
126     { "Spread Q Bspline",                2 }, 
127     { "Gather F",                      23  },
128     { "Gather F Bspline",              6   }, 
129     { "3D-FFT",                        8   },
130     { "Convolution",                   4   },
131     { "Solve PME",                     64  },
132     { "NS-Pairs",                      21  },
133     { "Reset In Box",                  3   },
134     { "Shift-X",                       6   },
135     { "CG-CoM",                        3   },
136     { "Sum Forces",                    1   },
137     { "Bonds",                         59  },
138     { "G96Bonds",                      44  },
139     { "FENE Bonds",                    58  },
140     { "Tab. Bonds",                    62  },
141     { "Restraint Potential",           86  },
142     { "Linear Angles",                 57  },
143     { "Angles",                        168 },
144     { "G96Angles",                     150 },
145     { "Quartic Angles",                160 },
146     { "Tab. Angles",                   169 },
147     { "Propers",                       229 },
148     { "Impropers",                     208 },
149     { "RB-Dihedrals",                  247 },
150     { "Four. Dihedrals",               247 },
151     { "Tab. Dihedrals",                227 },
152     { "Dist. Restr.",                  200 },
153     { "Orient. Restr.",                200 },
154     { "Dihedral Restr.",               200 },
155     { "Pos. Restr.",                   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
190 void init_nrnb(t_nrnb *nrnb)
191 {
192   int i;
193
194   for(i=0; (i<eNRNB); i++)
195     nrnb->n[i]=0.0;
196 }
197
198 void cp_nrnb(t_nrnb *dest, t_nrnb *src)
199 {
200   int i;
201
202   for(i=0; (i<eNRNB); i++)
203     dest->n[i]=src->n[i];
204 }
205
206 void add_nrnb(t_nrnb *dest, t_nrnb *s1, t_nrnb *s2)
207 {
208   int i;
209
210   for(i=0; (i<eNRNB); i++)
211     dest->n[i]=s1->n[i]+s2->n[i];
212 }
213
214 void print_nrnb(FILE *out, t_nrnb *nrnb)
215 {
216   int i;
217
218   for(i=0; (i<eNRNB); i++)
219     if (nrnb->n[i] > 0)
220       fprintf(out," %-26s %10.0f.\n",nbdata[i].name,nrnb->n[i]);
221 }
222
223 void _inc_nrnb(t_nrnb *nrnb,int enr,int inc,char *file,int line)
224 {
225   nrnb->n[enr]+=inc;
226 #ifdef DEBUG_NRNB
227   printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
228           nbdata[enr].name,enr,inc,file,line);
229 #endif
230 }
231
232 void print_flop(FILE *out,t_nrnb *nrnb,double *nbfs,double *mflop)
233 {
234   int    i;
235   double mni,frac,tfrac,tflop;
236   const char   *myline = "-----------------------------------------------------------------------------";
237   
238   *nbfs = 0.0;
239   for(i=0; (i<eNR_NBKERNEL_ALLVSALLGB); i++) {
240     if (strstr(nbdata[i].name,"W3-W3") != NULL)
241       *nbfs += 9e-6*nrnb->n[i];
242     else if (strstr(nbdata[i].name,"W3") != NULL)
243       *nbfs += 3e-6*nrnb->n[i];
244     else if (strstr(nbdata[i].name,"W4-W4") != NULL)
245       *nbfs += 10e-6*nrnb->n[i];
246     else if (strstr(nbdata[i].name,"W4") != NULL)
247       *nbfs += 4e-6*nrnb->n[i];
248     else
249       *nbfs += 1e-6*nrnb->n[i];
250   }
251   tflop=0;
252   for(i=0; (i<eNRNB); i++) 
253     tflop+=1e-6*nrnb->n[i]*nbdata[i].flop;
254   
255   if (tflop == 0) {
256     fprintf(out,"No MEGA Flopsen this time\n");
257     return;
258   }
259   if (out) {
260     fprintf(out,"\n\tM E G A - F L O P S   A C C O U N T I N G\n\n");
261   }
262
263   if (out)
264   {
265       fprintf(out," NB=Group-cutoff nonbonded kernels    NxN=N-by-N cluster Verlet kernels\n");
266       fprintf(out," RF=Reaction-Field  VdW=Van der Waals  QSTab=quadratic-spline table\n");
267       fprintf(out," W3=SPC/TIP3p  W4=TIP4p (single or pairs)\n");
268       fprintf(out," V&F=Potential and force  V=Potential only  F=Force only\n\n");
269
270       fprintf(out," %-32s %16s %15s  %7s\n",
271               "Computing:","M-Number","M-Flops","% Flops");
272       fprintf(out,"%s\n",myline);
273   }
274   *mflop=0.0;
275   tfrac=0.0;
276   for(i=0; (i<eNRNB); i++) {
277     mni     = 1e-6*nrnb->n[i];
278     *mflop += mni*nbdata[i].flop;
279     frac    = 100.0*mni*nbdata[i].flop/tflop;
280     tfrac  += frac;
281     if (out && mni != 0)
282       fprintf(out," %-32s %16.6f %15.3f  %6.1f\n",
283               nbdata[i].name,mni,mni*nbdata[i].flop,frac);
284   }
285   if (out) {
286     fprintf(out,"%s\n",myline);
287     fprintf(out," %-32s %16s %15.3f  %6.1f\n",
288             "Total","",*mflop,tfrac);
289     fprintf(out,"%s\n\n",myline);
290   }
291 }
292
293 void print_perf(FILE *out,double nodetime,double realtime,int nprocs,
294                 gmx_large_int_t nsteps,real delta_t,
295                 double nbfs,double mflop,
296                 int omp_nth_pp)
297 {
298   real runtime;
299
300   fprintf(out,"\n");
301
302   if (realtime > 0) 
303   {
304     fprintf(out,"%12s %12s %12s %10s\n","","Core t (s)","Wall t (s)","(%)");
305     fprintf(out,"%12s %12.3f %12.3f %10.1f\n","Time:",
306             nodetime, realtime, 100.0*nodetime/realtime);
307     /* only print day-hour-sec format if realtime is more than 30 min */
308     if (realtime > 30*60)
309     {
310       fprintf(out,"%12s %12s","","");
311       pr_difftime(out,realtime);
312     }
313     if (delta_t > 0) 
314     {
315       mflop = mflop/realtime;
316       runtime = nsteps*delta_t;
317
318       if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
319       {
320           fprintf(out,"%12s %12s %12s\n",
321                   "","(ns/day)","(hour/ns)");
322           fprintf(out,"%12s %12.3f %12.3f\n","Performance:",
323                   runtime*24*3.6/realtime,1000*realtime/(3600*runtime));
324       }
325       else
326       {
327         fprintf(out,"%12s %12s %12s %12s %12s\n",
328                 "","(Mnbf/s)",(mflop > 1000) ? "(GFlops)" : "(MFlops)",
329                 "(ns/day)","(hour/ns)");
330         fprintf(out,"%12s %12.3f %12.3f %12.3f %12.3f\n","Performance:",
331                 nbfs/realtime,(mflop > 1000) ? (mflop/1000) : mflop,
332                 runtime*24*3.6/realtime,1000*realtime/(3600*runtime));
333       }
334     } 
335     else 
336     {
337       if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
338       {
339           fprintf(out,"%12s %14s\n",
340                   "","(steps/hour)");
341           fprintf(out,"%12s %14.1f\n","Performance:",
342                   nsteps*3600.0/realtime);
343       }
344       else
345       {
346           fprintf(out,"%12s %12s %12s %14s\n",
347                   "","(Mnbf/s)",(mflop > 1000) ? "(GFlops)" : "(MFlops)",
348                   "(steps/hour)");
349           fprintf(out,"%12s %12.3f %12.3f %14.1f\n","Performance:",
350               nbfs/realtime,(mflop > 1000) ? (mflop/1000) : mflop,
351               nsteps*3600.0/realtime);
352       }
353     }
354   }
355 }
356
357 int cost_nrnb(int enr)
358 {
359   return nbdata[enr].flop;
360 }
361
362 const char *nrnb_str(int enr)
363 {
364   return nbdata[enr].name;
365 }
366
367 static const int    force_index[]={ 
368   eNR_BONDS,  eNR_ANGLES,  eNR_PROPER, eNR_IMPROPER, 
369   eNR_RB,     eNR_DISRES,  eNR_ORIRES, eNR_POSRES,
370   eNR_NS,
371 };
372 #define NFORCE_INDEX asize(force_index)
373
374 static const int    constr_index[]={ 
375   eNR_SHAKE,     eNR_SHAKE_RIJ, eNR_SETTLE,       eNR_UPDATE,       eNR_PCOUPL,
376   eNR_CONSTR_VIR,eNR_CONSTR_V
377 };
378 #define NCONSTR_INDEX asize(constr_index)
379
380 static double pr_av(FILE *log,t_commrec *cr,
381                     double fav,double ftot[],const char *title)
382 {
383   int    i,perc;
384   double dperc,unb;
385   
386   unb=0;
387   if (fav > 0) {
388     fav /= cr->nnodes - cr->npmenodes;
389     fprintf(log,"\n %-26s",title);
390     for(i=0; (i<cr->nnodes); i++) {
391         dperc=(100.0*ftot[i])/fav;
392         unb=max(unb,dperc);
393         perc=dperc;
394         fprintf(log,"%3d ",perc);
395     }
396     if (unb > 0) {
397       perc=10000.0/unb;
398       fprintf(log,"%6d%%\n\n",perc);
399     }
400     else
401       fprintf(log,"\n\n");
402   }
403   return unb;
404 }
405
406 void pr_load(FILE *log,t_commrec *cr,t_nrnb nrnb[])
407 {
408   int    i,j,perc;
409   double dperc,unb,uf,us;
410   double *ftot,fav;
411   double *stot,sav;
412   t_nrnb *av;
413
414   snew(av,1);
415   snew(ftot,cr->nnodes);
416   snew(stot,cr->nnodes);
417   init_nrnb(av);
418   for(i=0; (i<cr->nnodes); i++) {
419       add_nrnb(av,av,&(nrnb[i]));
420       /* Cost due to forces */
421       for(j=0; (j<eNR_NBKERNEL_ALLVSALLGB); j++)
422         ftot[i]+=nrnb[i].n[j]*cost_nrnb(j);
423       for(j=0; (j<NFORCE_INDEX); j++) 
424         ftot[i]+=nrnb[i].n[force_index[j]]*cost_nrnb(force_index[j]);
425       /* Due to shake */
426       for(j=0; (j<NCONSTR_INDEX); j++) {
427         stot[i]+=nrnb[i].n[constr_index[j]]*cost_nrnb(constr_index[j]);
428       }
429   }   
430   for(j=0; (j<eNRNB); j++)
431     av->n[j]=av->n[j]/(double)(cr->nnodes - cr->npmenodes);
432     
433     fprintf(log,"\nDetailed load balancing info in percentage of average\n");
434   
435   fprintf(log," Type                 NODE:");
436   for(i=0; (i<cr->nnodes); i++)
437       fprintf(log,"%3d ",i);
438   fprintf(log,"Scaling\n");
439   fprintf(log,"---------------------------");
440   for(i=0; (i<cr->nnodes); i++)
441       fprintf(log,"----");
442   fprintf(log,"-------\n");
443   
444   for(j=0; (j<eNRNB); j++) {
445     unb=100.0;
446     if (av->n[j] > 0) {
447       fprintf(log," %-26s",nrnb_str(j));
448       for(i=0; (i<cr->nnodes); i++) {
449           dperc=(100.0*nrnb[i].n[j])/av->n[j];
450           unb=max(unb,dperc);
451           perc=dperc;
452           fprintf(log,"%3d ",perc);
453       }
454       if (unb > 0) {
455         perc=10000.0/unb;
456         fprintf(log,"%6d%%\n",perc);
457       }
458       else
459         fprintf(log,"\n");
460     }   
461   }
462   fav=sav=0;
463   for(i=0; (i<cr->nnodes); i++) {
464     fav+=ftot[i];
465     sav+=stot[i];
466   }
467   uf=pr_av(log,cr,fav,ftot,"Total Force");
468   us=pr_av(log,cr,sav,stot,"Total Constr.");
469   
470   unb=(uf*fav+us*sav)/(fav+sav);
471   if (unb > 0) {
472     unb=10000.0/unb;
473     fprintf(log,"\nTotal Scaling: %.0f%% of max performance\n\n",unb);
474   }
475 }
476