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