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