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