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