Merge release-4-6 into master
[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 #include "copyrite.h"
49
50 typedef struct {
51   const char *name;
52   int  flop;
53 } t_nrnb_data;
54
55
56 static const t_nrnb_data nbdata[eNRNB] = {
57     { "LJ",                             33 }, /* nb_kernel010 */
58     { "Buckingham",                     61 }, /* nb_kernel020 */ 
59     { "VdW(T)",                         54 }, /* nb_kernel030 */
60     { "Coulomb",                        27 }, /* nb_kernel100 */
61     { "Coulomb [W3]",                   80 }, /* nb_kernel101 */
62     { "Coulomb [W3-W3]",               234 }, /* nb_kernel102 */
63     { "Coulomb [W4]",                   80 }, /* nb_kernel103 */
64     { "Coulomb [W4-W4]",               234 }, /* nb_kernel104 */
65     { "Coulomb + LJ",                   38 }, /* nb_kernel110 */
66     { "Coulomb + LJ [W3]",              91 }, /* nb_kernel111 */
67     { "Coulomb + LJ [W3-W3]",          245 }, /* nb_kernel112 */
68     { "Coulomb + LJ [W4]",             113 }, /* nb_kernel113 */
69     { "Coulomb + LJ [W4-W4]",          267 }, /* nb_kernel114 */
70     { "Coulomb + Bham ",                64 }, /* nb_kernel120 */
71     { "Coulomb + Bham [W3]",           117 }, /* nb_kernel121 */
72     { "Coulomb + Bham [W3-W3]",        271 }, /* nb_kernel122 */
73     { "Coulomb + Bham [W4]",           141 }, /* nb_kernel123 */
74     { "Coulomb + Bham [W4-W4]",        295 }, /* nb_kernel124 */
75     { "Coulomb + VdW(T) ",              59 }, /* nb_kernel130 */
76     { "Coulomb + VdW(T) [W3]",         112 }, /* nb_kernel131 */
77     { "Coulomb + VdW(T) [W3-W3]",      266 }, /* nb_kernel132 */
78     { "Coulomb + VdW(T) [W4]",         134 }, /* nb_kernel133 */
79     { "Coulomb + VdW(T) [W4-W4]",      288 }, /* nb_kernel134 */
80     { "RF Coul",                        33 }, /* nb_kernel200 */
81     { "RF Coul [W3]",                   98 }, /* nb_kernel201 */
82     { "RF Coul [W3-W3]",               288 }, /* nb_kernel202 */
83     { "RF Coul [W4]",                   98 }, /* nb_kernel203 */
84     { "RF Coul [W4-W4]",               288 }, /* nb_kernel204 */
85     { "RF Coul + LJ",                   44 }, /* nb_kernel210 */
86     { "RF Coul + LJ [W3]",             109 }, /* nb_kernel211 */
87     { "RF Coul + LJ [W3-W3]",          299 }, /* nb_kernel212 */
88     { "RF Coul + LJ [W4]",             131 }, /* nb_kernel213 */
89     { "RF Coul + LJ [W4-W4]",          321 }, /* nb_kernel214 */
90     { "RF Coul + Bham ",                70 }, /* nb_kernel220 */
91     { "RF Coul + Bham [W3]",           135 }, /* nb_kernel221 */
92     { "RF Coul + Bham [W3-W3]",        325 }, /* nb_kernel222 */
93     { "RF Coul + Bham [W4]",           159 }, /* nb_kernel223 */
94     { "RF Coul + Bham [W4-W4]",        349 }, /* nb_kernel224 */
95     { "RF Coul + VdW(T) ",              65 }, /* nb_kernel230 */
96     { "RF Coul + VdW(T) [W3]",         130 }, /* nb_kernel231 */
97     { "RF Coul + VdW(T) [W3-W3]",      320 }, /* nb_kernel232 */
98     { "RF Coul + VdW(T) [W4]",         152 }, /* nb_kernel233 */
99     { "RF Coul + VdW(T) [W4-W4]",      342 }, /* nb_kernel234 */
100     { "Coul(T)",                        42 }, /* nb_kernel300 */
101     { "Coul(T) [W3]",                  125 }, /* nb_kernel301 */
102     { "Coul(T) [W3-W3]",               369 }, /* nb_kernel302 */
103     { "Coul(T) [W4]",                  125 }, /* nb_kernel303 */
104     { "Coul(T) [W4-W4]",               369 }, /* nb_kernel304 */
105     { "Coul(T) + LJ",                   55 }, /* nb_kernel310 */
106     { "Coul(T) + LJ [W3]",             138 }, /* nb_kernel311 */
107     { "Coul(T) + LJ [W3-W3]",          382 }, /* nb_kernel312 */
108     { "Coul(T) + LJ [W4]",             158 }, /* nb_kernel313 */
109     { "Coul(T) + LJ [W4-W4]",          402 }, /* nb_kernel314 */
110     { "Coul(T) + Bham",                 81 }, /* nb_kernel320 */
111     { "Coul(T) + Bham [W3]",           164 }, /* nb_kernel321 */
112     { "Coul(T) + Bham [W3-W3]",        408 }, /* nb_kernel322 */
113     { "Coul(T) + Bham [W4]",           186 }, /* nb_kernel323 */
114     { "Coul(T) + Bham [W4-W4]",        430 }, /* nb_kernel324 */
115     { "Coul(T) + VdW(T)",               68 }, /* nb_kernel330 */
116     { "Coul(T) + VdW(T) [W3]",         151 }, /* nb_kernel331 */
117     { "Coul(T) + VdW(T) [W3-W3]",      395 }, /* nb_kernel332 */
118     { "Coul(T) + VdW(T) [W4]",         179 }, /* nb_kernel333 */
119     { "Coul(T) + VdW(T) [W4-W4]",      423 }, /* nb_kernel334 */
120     { "Generalized Born Coulomb",       48 }, /* nb_kernel400 */
121     { "GB Coulomb + LJ",                61 }, /* nb_kernel410 */
122     { "GB Coulomb + VdW(T)",            79 }, /* nb_kernel430 */
123     { "LJ NF",                          19 }, /* nb_kernel010nf */
124     { "Buckingham NF",                  48 }, /* nb_kernel020nf */ 
125     { "VdW(T) NF",                      33 }, /* nb_kernel030nf */
126     { "Coulomb NF",                     16 }, /* nb_kernel100nf */
127     { "Coulomb [W3] NF",                47 }, /* nb_kernel101nf */
128     { "Coulomb [W3-W3] NF",            135 }, /* nb_kernel102nf */
129     { "Coulomb [W4] NF",                47 }, /* nb_kernel103nf */
130     { "Coulomb [W4-W4] NF",            135 }, /* nb_kernel104nf */
131     { "Coulomb + LJ NF",                24 }, /* nb_kernel110nf */
132     { "Coulomb + LJ [W3] NF",           55 }, /* nb_kernel111nf */
133     { "Coulomb + LJ [W3-W3] NF",       143 }, /* nb_kernel112nf */
134     { "Coulomb + LJ [W4] NF",           66 }, /* nb_kernel113nf */
135     { "Coulomb + LJ [W4-W4] NF",       154 }, /* nb_kernel114nf */
136     { "Coulomb + Bham  NF",             51 }, /* nb_kernel120nf */
137     { "Coulomb + Bham [W3] NF",         82 }, /* nb_kernel121nf */
138     { "Coulomb + Bham [W3-W3] NF",     170 }, /* nb_kernel122nf */
139     { "Coulomb + Bham [W4] NF",         95 }, /* nb_kernel123nf */
140     { "Coulomb + Bham [W4-W4] NF",     183 }, /* nb_kernel124nf */
141     { "Coulomb + VdW(T)  NF",           36 }, /* nb_kernel130nf */
142     { "Coulomb + VdW(T) [W3] NF",       67 }, /* nb_kernel131nf */
143     { "Coulomb + VdW(T) [W3-W3] NF",   155 }, /* nb_kernel132nf */
144     { "Coulomb + VdW(T) [W4] NF",       80 }, /* nb_kernel133nf */
145     { "Coulomb + VdW(T) [W4-W4] NF",   168 }, /* nb_kernel134nf */
146     { "RF Coul NF",                     19 }, /* nb_kernel200nf */
147     { "RF Coul [W3] NF",                56 }, /* nb_kernel201nf */
148     { "RF Coul [W3-W3] NF",            162 }, /* nb_kernel202nf */
149     { "RF Coul [W4] NF",                56 }, /* nb_kernel203nf */
150     { "RF Coul [W4-W4] NF",            162 }, /* nb_kernel204nf */
151     { "RF Coul + LJ NF",                27 }, /* nb_kernel210nf */
152     { "RF Coul + LJ [W3] NF",           64 }, /* nb_kernel211nf */
153     { "RF Coul + LJ [W3-W3] NF",       170 }, /* nb_kernel212nf */
154     { "RF Coul + LJ [W4] NF",           75 }, /* nb_kernel213nf */
155     { "RF Coul + LJ [W4-W4] NF",       181 }, /* nb_kernel214nf */
156     { "RF Coul + Bham  NF",             54 }, /* nb_kernel220nf */
157     { "RF Coul + Bham [W3] NF",         91 }, /* nb_kernel221nf */
158     { "RF Coul + Bham [W3-W3] NF",     197 }, /* nb_kernel222nf */
159     { "RF Coul + Bham [W4] NF",        104 }, /* nb_kernel223nf */
160     { "RF Coul + Bham [W4-W4] NF",     210 }, /* nb_kernel224nf */
161     { "RF Coul + VdW(T)  NF",           39 }, /* nb_kernel230nf */
162     { "RF Coul + VdW(T) [W3] NF",       76 }, /* nb_kernel231nf */
163     { "RF Coul + VdW(T) [W3-W3] NF",   182 }, /* nb_kernel232nf */
164     { "RF Coul + VdW(T) [W4] NF",       89 }, /* nb_kernel233nf */
165     { "RF Coul + VdW(T) [W4-W4] NF",   195 }, /* nb_kernel234nf */
166     { "Coul(T) NF",                     26 }, /* nb_kernel300nf */
167     { "Coul(T) [W3] NF",                77 }, /* nb_kernel301nf */
168     { "Coul(T) [W3-W3] NF",            225 }, /* nb_kernel302nf */
169     { "Coul(T) [W4] NF",                77 }, /* nb_kernel303nf */
170     { "Coul(T) [W4-W4] NF",            225 }, /* nb_kernel304nf */
171     { "Coul(T) + LJ NF",                34 }, /* nb_kernel310nf */
172     { "Coul(T) + LJ [W3] NF",           85 }, /* nb_kernel311nf */
173     { "Coul(T) + LJ [W3-W3] NF",       233 }, /* nb_kernel312nf */
174     { "Coul(T) + LJ [W4] NF",           96 }, /* nb_kernel313nf */
175     { "Coul(T) + LJ [W4-W4] NF",       244 }, /* nb_kernel314nf */
176     { "Coul(T) + Bham NF",              61 }, /* nb_kernel320nf */
177     { "Coul(T) + Bham [W3] NF",        112 }, /* nb_kernel321nf */
178     { "Coul(T) + Bham [W3-W3] NF",     260 }, /* nb_kernel322nf */
179     { "Coul(T) + Bham [W4] NF",        125 }, /* nb_kernel323nf */
180     { "Coul(T) + Bham [W4-W4] NF",     273 }, /* nb_kernel324nf */
181     { "Coul(T) + VdW(T) NF",            42 }, /* nb_kernel330nf */
182     { "Coul(T) + VdW(T) [W3] NF",       93 }, /* nb_kernel331nf */
183     { "Coul(T) + VdW(T) [W3-W3] NF",   241 }, /* nb_kernel332nf */
184     { "Coul(T) + VdW(T) [W4] NF",      110 }, /* nb_kernel333nf */
185     { "Coul(T) + VdW(T) [W4-W4] NF",   258 }, /* nb_kernel334nf */
186     { "Generalized Born Coulomb NF",    29 }, /* nb_kernel400nf */
187     { "GB Coulomb + LJ NF",             37 }, /* nb_kernel410nf */
188     { "GB Coulomb + VdW(T) NF",         49 }, /* nb_kernel430nf */
189     { "Free energy innerloop",         150 }, /* free energy, estimate */  
190     { "All-vs-All, Coul + LJ",          38 },
191     { "All-vs-All, GB + LJ",            61 },
192     { "Outer nonbonded loop",           10 },
193     { "Pair Search distance check",      9 }, /* nbnxn pair dist. check */
194     /* nbnxn kernel flops are based on inner-loops without exclusion checks.
195      * Plain Coulomb runs through the RF kernels, except with CUDA.
196      * invsqrt is counted as 6 flops: 1 for _mm_rsqt_ps + 5 for iteration.
197      * The flops are equal for plain-C, x86 SIMD and CUDA, except for:
198      * - plain-C kernel uses one flop more for Coulomb-only (F) than listed
199      * - x86 SIMD LJ geom-comb.rule kernels (fastest) use 2 more flops
200      * - x86 SIMD LJ LB-comb.rule kernels (fast) use 3 (8 for F+E) more flops
201      * - GPU always does exclusions, which requires 2-4 flops, but as invsqrt
202      *   is always counted as 6 flops, this roughly compensates.
203      */
204     { "LJ + Coulomb RF (F)",            38 }, /* nbnxn kernel LJ+RF, no ener */
205     { "LJ + Coulomb RF (F+E)",          54 },
206     { "LJ + Coulomb tabulated (F)",     41 }, /* nbnxn kernel LJ+tab, no en */
207     { "LJ + Coulomb tabulated (F+E)",   59 },
208     { "LJ (F)",                         33 }, /* nbnxn kernel LJ, no ener */
209     { "LJ (F+E)",                       43 },
210     { "Coulomb RF (F)",                 31 }, /* nbnxn kernel RF, no ener */
211     { "Coulomb RF (F+E)",               36 },
212     { "Coulomb tabulated (F)",          34 }, /* nbnxn kernel tab, no ener */
213     { "Coulomb tabulated (F+E)",        41 },
214     { "1,4 nonbonded interactions",     90 },
215     { "Born radii (Still)",             47 },
216     { "Born radii (HCT/OBC)",          183 },
217     { "Born force chain rule",          15 },
218     { "All-vs-All Still radii",         47 },
219     { "All-vs-All HCT/OBC radii",      183 },
220     { "All-vs-All Born chain rule",     15 },
221     { "Calc Weights",                   36 },
222     { "Spread Q",                        6 },
223     { "Spread Q Bspline",                2 }, 
224     { "Gather F",                      23  },
225     { "Gather F Bspline",              6   }, 
226     { "3D-FFT",                        8   },
227     { "Convolution",                   4   },
228     { "Solve PME",                     64  },
229     { "NS-Pairs",                      21  },
230     { "Reset In Box",                  3   },
231     { "Shift-X",                       6   },
232     { "CG-CoM",                        3   },
233     { "Sum Forces",                    1   },
234     { "Bonds",                         59  },
235     { "G96Bonds",                      44  },
236     { "FENE Bonds",                    58  },
237     { "Tab. Bonds",                    62  },
238     { "Restraint Potential",           86  },
239     { "Linear Angles",                 57  },
240     { "Angles",                        168 },
241     { "G96Angles",                     150 },
242     { "Quartic Angles",                160 },
243     { "Tab. Angles",                   169 },
244     { "Propers",                       229 },
245     { "Impropers",                     208 },
246     { "RB-Dihedrals",                  247 },
247     { "Four. Dihedrals",               247 },
248     { "Tab. Dihedrals",                227 },
249     { "Dist. Restr.",                  200 },
250     { "Orient. Restr.",                200 },
251     { "Dihedral Restr.",               200 },
252     { "Pos. Restr.",                   50  },
253     { "Flat-bottom posres",            50  },
254     { "Angle Restr.",                  191 },
255     { "Angle Restr. Z",                164 },
256     { "Morse Potent.",                 83  },
257     { "Cubic Bonds",                   54  },
258     { "Walls",                         31  },
259     { "Polarization",                  59  },
260     { "Anharmonic Polarization",       72  },
261     { "Water Pol.",                    62  },
262     { "Thole Pol.",                    296 },
263     { "Virial",                        18  },
264     { "Update",                        31  },
265     { "Ext.ens. Update",               54  },
266     { "Stop-CM",                       10  },
267     { "P-Coupling",                    6   },
268     { "Calc-Ekin",                     27  },
269     { "Lincs",                         60  },
270     { "Lincs-Mat",                     4   },
271     { "Shake",                         30  },
272     { "Constraint-V",                   8  },
273     { "Shake-Init",                    10  },
274     { "Constraint-Vir",                24  },
275     { "Settle",                        323 },
276     { "Virtual Site 2",                23  },
277     { "Virtual Site 3",                37  },
278     { "Virtual Site 3fd",              95  },
279     { "Virtual Site 3fad",             176 },
280     { "Virtual Site 3out",             87  },
281     { "Virtual Site 4fd",              110 }, 
282     { "Virtual Site 4fdn",             254 }, 
283     { "Virtual Site N",                 15 },
284     { "Mixed Generalized Born stuff",   10 } 
285 };
286
287
288 void init_nrnb(t_nrnb *nrnb)
289 {
290   int i;
291
292   for(i=0; (i<eNRNB); i++)
293     nrnb->n[i]=0.0;
294 }
295
296 void cp_nrnb(t_nrnb *dest, t_nrnb *src)
297 {
298   int i;
299
300   for(i=0; (i<eNRNB); i++)
301     dest->n[i]=src->n[i];
302 }
303
304 void add_nrnb(t_nrnb *dest, t_nrnb *s1, t_nrnb *s2)
305 {
306   int i;
307
308   for(i=0; (i<eNRNB); i++)
309     dest->n[i]=s1->n[i]+s2->n[i];
310 }
311
312 void print_nrnb(FILE *out, t_nrnb *nrnb)
313 {
314   int i;
315
316   for(i=0; (i<eNRNB); i++)
317     if (nrnb->n[i] > 0)
318       fprintf(out," %-26s %10.0f.\n",nbdata[i].name,nrnb->n[i]);
319 }
320
321 void _inc_nrnb(t_nrnb *nrnb,int enr,int inc,char *file,int line)
322 {
323   nrnb->n[enr]+=inc;
324 #ifdef DEBUG_NRNB
325   printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
326           nbdata[enr].name,enr,inc,file,line);
327 #endif
328 }
329
330 void print_flop(FILE *out,t_nrnb *nrnb,double *nbfs,double *mflop)
331 {
332   int    i;
333   double mni,frac,tfrac,tflop;
334   const char   *myline = "-----------------------------------------------------------------------------";
335   
336   *nbfs = 0.0;
337   for(i=0; (i<eNR_NBKERNEL_NR); i++) {
338     if (strstr(nbdata[i].name,"W3-W3") != NULL)
339       *nbfs += 9e-6*nrnb->n[i];
340     else if (strstr(nbdata[i].name,"W3") != NULL)
341       *nbfs += 3e-6*nrnb->n[i];
342     else if (strstr(nbdata[i].name,"W4-W4") != NULL)
343       *nbfs += 10e-6*nrnb->n[i];
344     else if (strstr(nbdata[i].name,"W4") != NULL)
345       *nbfs += 4e-6*nrnb->n[i];
346     else
347       *nbfs += 1e-6*nrnb->n[i];
348   }
349   tflop=0;
350   for(i=0; (i<eNRNB); i++) 
351     tflop+=1e-6*nrnb->n[i]*nbdata[i].flop;
352   
353   if (tflop == 0) {
354     fprintf(out,"No MEGA Flopsen this time\n");
355     return;
356   }
357   if (out) {
358     fprintf(out,"\n\tM E G A - F L O P S   A C C O U N T I N G\n\n");
359   }
360
361   if (out) {
362     fprintf(out,"   RF=Reaction-Field  FE=Free Energy  SCFE=Soft-Core/Free Energy\n");
363     fprintf(out,"   T=Tabulated        W3=SPC/TIP3p    W4=TIP4p (single or pairs)\n");
364     fprintf(out,"   NF=No Forces\n\n");
365     
366     fprintf(out," %-32s %16s %15s  %7s\n",
367             "Computing:","M-Number","M-Flops","% Flops");
368     fprintf(out,"%s\n",myline);
369   }
370   *mflop=0.0;
371   tfrac=0.0;
372   for(i=0; (i<eNRNB); i++) {
373     mni     = 1e-6*nrnb->n[i];
374     *mflop += mni*nbdata[i].flop;
375     frac    = 100.0*mni*nbdata[i].flop/tflop;
376     tfrac  += frac;
377     if (out && mni != 0)
378       fprintf(out," %-32s %16.6f %15.3f  %6.1f\n",
379               nbdata[i].name,mni,mni*nbdata[i].flop,frac);
380   }
381   if (out) {
382     fprintf(out,"%s\n",myline);
383     fprintf(out," %-32s %16s %15.3f  %6.1f\n",
384             "Total","",*mflop,tfrac);
385     fprintf(out,"%s\n\n",myline);
386   }
387 }
388
389 void print_perf(FILE *out,double nodetime,double realtime,int nprocs,
390                 gmx_large_int_t nsteps,real delta_t,
391                 double nbfs,double mflop,
392                 int omp_nth_pp)
393 {
394   real runtime;
395
396   fprintf(out,"\n");
397
398   if (realtime > 0) 
399   {
400     fprintf(out,"%12s %12s %12s %10s\n","","Core t (s)","Wall t (s)","(%)");
401     fprintf(out,"%12s %12.3f %12.3f %10.1f\n","Time:",
402             nodetime, realtime, 100.0*nodetime/realtime);
403     /* only print day-hour-sec format if realtime is more than 30 min */
404     if (realtime > 30*60)
405     {
406       fprintf(out,"%12s %12s","","");
407       pr_difftime(out,realtime);
408     }
409     if (delta_t > 0) 
410     {
411       mflop = mflop/realtime;
412       runtime = nsteps*delta_t;
413
414       if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
415       {
416           fprintf(out,"%12s %12s %12s\n",
417                   "","(ns/day)","(hour/ns)");
418           fprintf(out,"%12s %12.3f %12.3f\n","Performance:",
419                   runtime*24*3.6/realtime,1000*realtime/(3600*runtime));
420       }
421       else
422       {
423         fprintf(out,"%12s %12s %12s %12s %12s\n",
424                 "","(Mnbf/s)",(mflop > 1000) ? "(GFlops)" : "(MFlops)",
425                 "(ns/day)","(hour/ns)");
426         fprintf(out,"%12s %12.3f %12.3f %12.3f %12.3f\n","Performance:",
427                 nbfs/realtime,(mflop > 1000) ? (mflop/1000) : mflop,
428                 runtime*24*3.6/realtime,1000*realtime/(3600*runtime));
429       }
430     } 
431     else 
432     {
433       if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
434       {
435           fprintf(out,"%12s %14s\n",
436                   "","(steps/hour)");
437           fprintf(out,"%12s %14.1f\n","Performance:",
438                   nsteps*3600.0/realtime);
439       }
440       else
441       {
442           fprintf(out,"%12s %12s %12s %14s\n",
443                   "","(Mnbf/s)",(mflop > 1000) ? "(GFlops)" : "(MFlops)",
444                   "(steps/hour)");
445           fprintf(out,"%12s %12.3f %12.3f %14.1f\n","Performance:",
446               nbfs/realtime,(mflop > 1000) ? (mflop/1000) : mflop,
447               nsteps*3600.0/realtime);
448       }
449     }
450   }
451 }
452
453 int cost_nrnb(int enr)
454 {
455   return nbdata[enr].flop;
456 }
457
458 const char *nrnb_str(int enr)
459 {
460   return nbdata[enr].name;
461 }
462
463 static const int    force_index[]={ 
464   eNR_BONDS,  eNR_ANGLES,  eNR_PROPER, eNR_IMPROPER, 
465   eNR_RB,     eNR_DISRES,  eNR_ORIRES, eNR_POSRES,
466   eNR_FBPOSRES,  eNR_NS,     eNR_NBKERNEL_OUTER
467 };
468 #define NFORCE_INDEX asize(force_index)
469
470 static const int    constr_index[]={ 
471   eNR_SHAKE,     eNR_SHAKE_RIJ, eNR_SETTLE,       eNR_UPDATE,       eNR_PCOUPL,
472   eNR_CONSTR_VIR,eNR_CONSTR_V
473 };
474 #define NCONSTR_INDEX asize(constr_index)
475
476 static double pr_av(FILE *log,t_commrec *cr,
477                     double fav,double ftot[],const char *title)
478 {
479   int    i,perc;
480   double dperc,unb;
481   
482   unb=0;
483   if (fav > 0) {
484     fav /= cr->nnodes - cr->npmenodes;
485     fprintf(log,"\n %-26s",title);
486     for(i=0; (i<cr->nnodes); i++) {
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       perc=10000.0/unb;
494       fprintf(log,"%6d%%\n\n",perc);
495     }
496     else
497       fprintf(log,"\n\n");
498   }
499   return unb;
500 }
501
502 void pr_load(FILE *log,t_commrec *cr,t_nrnb nrnb[])
503 {
504   int    i,j,perc;
505   double dperc,unb,uf,us;
506   double *ftot,fav;
507   double *stot,sav;
508   t_nrnb *av;
509
510   snew(av,1);
511   snew(ftot,cr->nnodes);
512   snew(stot,cr->nnodes);
513   init_nrnb(av);
514   for(i=0; (i<cr->nnodes); i++) {
515       add_nrnb(av,av,&(nrnb[i]));
516       /* Cost due to forces */
517       for(j=0; (j<eNR_NBKERNEL_NR); j++)
518         ftot[i]+=nrnb[i].n[j]*cost_nrnb(j);
519       for(j=0; (j<NFORCE_INDEX); j++) 
520         ftot[i]+=nrnb[i].n[force_index[j]]*cost_nrnb(force_index[j]);
521       /* Due to shake */
522       for(j=0; (j<NCONSTR_INDEX); j++) {
523         stot[i]+=nrnb[i].n[constr_index[j]]*cost_nrnb(constr_index[j]);
524       }
525   }   
526   for(j=0; (j<eNRNB); j++)
527     av->n[j]=av->n[j]/(double)(cr->nnodes - cr->npmenodes);
528     
529     fprintf(log,"\nDetailed load balancing info in percentage of average\n");
530   
531   fprintf(log," Type                 NODE:");
532   for(i=0; (i<cr->nnodes); i++)
533       fprintf(log,"%3d ",i);
534   fprintf(log,"Scaling\n");
535   fprintf(log,"---------------------------");
536   for(i=0; (i<cr->nnodes); i++)
537       fprintf(log,"----");
538   fprintf(log,"-------\n");
539   
540   for(j=0; (j<eNRNB); j++) {
541     unb=100.0;
542     if (av->n[j] > 0) {
543       fprintf(log," %-26s",nrnb_str(j));
544       for(i=0; (i<cr->nnodes); i++) {
545           dperc=(100.0*nrnb[i].n[j])/av->n[j];
546           unb=max(unb,dperc);
547           perc=dperc;
548           fprintf(log,"%3d ",perc);
549       }
550       if (unb > 0) {
551         perc=10000.0/unb;
552         fprintf(log,"%6d%%\n",perc);
553       }
554       else
555         fprintf(log,"\n");
556     }   
557   }
558   fav=sav=0;
559   for(i=0; (i<cr->nnodes); i++) {
560     fav+=ftot[i];
561     sav+=stot[i];
562   }
563   uf=pr_av(log,cr,fav,ftot,"Total Force");
564   us=pr_av(log,cr,sav,stot,"Total Constr.");
565   
566   unb=(uf*fav+us*sav)/(fav+sav);
567   if (unb > 0) {
568     unb=10000.0/unb;
569     fprintf(log,"\nTotal Scaling: %.0f%% of max performance\n\n",unb);
570   }
571 }
572