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