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