Split lines with many copyright years
[alexxy/gromacs.git] / src / gromacs / ewald / pme_load_balancing.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  *
38  * \brief This file contains function definitions necessary for
39  * managing automatic load balance of PME calculations (Coulomb and
40  * LJ).
41  *
42  * \author Berk Hess <hess@kth.se>
43  * \ingroup module_ewald
44  */
45 #include "gmxpre.h"
46
47 #include "pme_load_balancing.h"
48
49 #include <cassert>
50 #include <cmath>
51
52 #include <algorithm>
53
54 #include "gromacs/domdec/dlb.h"
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_network.h"
57 #include "gromacs/domdec/domdec_struct.h"
58 #include "gromacs/domdec/partition.h"
59 #include "gromacs/ewald/ewald_utils.h"
60 #include "gromacs/ewald/pme.h"
61 #include "gromacs/fft/calcgrid.h"
62 #include "gromacs/gmxlib/network.h"
63 #include "gromacs/math/functions.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/mdlib/dispersioncorrection.h"
66 #include "gromacs/mdlib/forcerec.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/nbnxm/gpu_data_mgmt.h"
72 #include "gromacs/nbnxm/nbnxm.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/logger.h"
79 #include "gromacs/utility/smalloc.h"
80 #include "gromacs/utility/strconvert.h"
81
82 #include "pme_internal.h"
83
84 /*! \brief Parameters and settings for one PP-PME setup */
85 struct pme_setup_t
86 {
87     real rcut_coulomb;         /**< Coulomb cut-off                              */
88     real rlistOuter;           /**< cut-off for the outer pair-list              */
89     real rlistInner;           /**< cut-off for the inner pair-list              */
90     real spacing;              /**< (largest) PME grid spacing                   */
91     ivec grid;                 /**< the PME grid dimensions                      */
92     real grid_efficiency;      /**< ineffiency factor for non-uniform grids <= 1 */
93     real ewaldcoeff_q;         /**< Electrostatic Ewald coefficient            */
94     real ewaldcoeff_lj;        /**< LJ Ewald coefficient, only for the call to send_switchgrid */
95     struct gmx_pme_t* pmedata; /**< the data structure used in the PME code      */
96     int               count;   /**< number of times this setup has been timed    */
97     double            cycles;  /**< the fastest time for this setup in cycles    */
98 };
99
100 /*! \brief After 50 nstlist periods of not observing imbalance: never tune PME */
101 const int PMETunePeriod = 50;
102 /*! \brief Trigger PME load balancing at more than 5% PME overload */
103 const real loadBalanceTriggerFactor = 1.05;
104 /*! \brief Scale the grid by a most at factor 1.7.
105  *
106  * This still leaves room for about 4-4.5x decrease in grid spacing while limiting the cases where
107  * large imbalance leads to extreme cutoff scaling for marginal benefits.
108  *
109  * This should help to avoid:
110  *   - large increase in power consumption for little performance gain
111  *   - increasing communication volume
112  *   - limiting DLB
113  */
114 const real c_maxSpacingScaling = 1.7;
115 /*! \brief In the initial scan, step by grids that are at least a factor 0.8 coarser */
116 const real gridpointsScaleFactor = 0.8;
117 /*! \brief In the initial scan, try to skip grids with uneven x/y/z spacing,
118  * checking if the "efficiency" is more than 5% worse than the previous grid.
119  */
120 const real relativeEfficiencyFactor = 1.05;
121 /*! \brief Rerun until a run is 12% slower setups than the fastest run so far */
122 const real maxRelativeSlowdownAccepted = 1.12;
123 /*! \brief If setups get more than 2% faster, do another round to avoid
124  * choosing a slower setup due to acceleration or fluctuations.
125  */
126 const real maxFluctuationAccepted = 1.02;
127
128 //! \brief Number of nstlist long tuning intervals to skip before starting
129 //         load-balancing at the beginning of the run.
130 const int c_numFirstTuningIntervalSkip = 5;
131 //! \brief Number of nstlist long tuning intervals to skip before starting
132 //         load-balancing at the beginning of the run with separate PME ranks. */
133 const int c_numFirstTuningIntervalSkipWithSepPme = 3;
134 //! \brief Number of nstlist long tuning intervals to skip after switching to a new setting
135 //         during balancing.
136 const int c_numPostSwitchTuningIntervalSkip = 1;
137 //! \brief Number of seconds to delay the tuning at startup to allow processors clocks to ramp up.
138 const double c_startupTimeDelay = 5.0;
139
140 /*! \brief Enumeration whose values describe the effect limiting the load balancing */
141 enum epmelb
142 {
143     epmelblimNO,
144     epmelblimBOX,
145     epmelblimDD,
146     epmelblimPMEGRID,
147     epmelblimMAXSCALING,
148     epmelblimNR
149 };
150
151 /*! \brief Descriptive strings matching ::epmelb */
152 static const char* pmelblim_str[epmelblimNR] = { "no", "box size", "domain decompostion",
153                                                  "PME grid restriction",
154                                                  "maximum allowed grid scaling" };
155
156 struct pme_load_balancing_t
157 {
158     gmx_bool bSepPMERanks;  /**< do we have separate PME ranks? */
159     gmx_bool bActive;       /**< is PME tuning active? */
160     int64_t  step_rel_stop; /**< stop the tuning after this value of step_rel */
161     gmx_bool bTriggerOnDLB; /**< trigger balancing only on DD DLB */
162     gmx_bool bBalance;      /**< are we in the balancing phase, i.e. trying different setups? */
163     int      nstage;        /**< the current maximum number of stages */
164
165     real                     cut_spacing;        /**< the minimum cutoff / PME grid spacing ratio */
166     real                     rcut_vdw;           /**< Vdw cutoff (does not change) */
167     real                     rcut_coulomb_start; /**< Initial electrostatics cutoff */
168     real                     rbufOuter_coulomb;  /**< the outer pairlist buffer size */
169     real                     rbufOuter_vdw;      /**< the outer pairlist buffer size */
170     real                     rbufInner_coulomb;  /**< the inner pairlist buffer size */
171     real                     rbufInner_vdw;      /**< the inner pairlist buffer size */
172     matrix                   box_start;          /**< the initial simulation box */
173     std::vector<pme_setup_t> setup;              /**< the PME+cutoff setups */
174     int                      cur;                /**< the index (in setup) of the current setup */
175     int                      fastest;            /**< index of the fastest setup up till now */
176     int                      lower_limit;        /**< don't go below this setup index */
177     int                      start;    /**< start of setup index range to consider in stage>0 */
178     int                      end;      /**< end   of setup index range to consider in stage>0 */
179     int                      elimited; /**< was the balancing limited, uses enum above */
180     int                      cutoff_scheme; /**< Verlet or group cut-offs */
181
182     int stage; /**< the current stage */
183
184     int    cycles_n; /**< step cycle counter cumulative count */
185     double cycles_c; /**< step cycle counter cumulative cycles */
186     double startTime; /**< time stamp when the balancing was started (relative to the UNIX epoch start).*/
187 };
188
189 /* TODO The code in this file should call this getter, rather than
190  * read bActive anywhere */
191 bool pme_loadbal_is_active(const pme_load_balancing_t* pme_lb)
192 {
193     return pme_lb != nullptr && pme_lb->bActive;
194 }
195
196 // TODO Return a unique_ptr to pme_load_balancing_t
197 void pme_loadbal_init(pme_load_balancing_t**     pme_lb_p,
198                       t_commrec*                 cr,
199                       const gmx::MDLogger&       mdlog,
200                       const t_inputrec&          ir,
201                       const matrix               box,
202                       const interaction_const_t& ic,
203                       const nonbonded_verlet_t&  nbv,
204                       gmx_pme_t*                 pmedata,
205                       gmx_bool                   bUseGPU)
206 {
207
208     pme_load_balancing_t* pme_lb;
209     real                  spm, sp;
210     int                   d;
211
212     // Note that we don't (yet) support PME load balancing with LJ-PME only.
213     GMX_RELEASE_ASSERT(EEL_PME(ir.coulombtype),
214                        "pme_loadbal_init called without PME electrostatics");
215     // To avoid complexity, we require a single cut-off with PME for q+LJ.
216     // This is checked by grompp, but it doesn't hurt to check again.
217     GMX_RELEASE_ASSERT(!(EEL_PME(ir.coulombtype) && EVDW_PME(ir.vdwtype) && ir.rcoulomb != ir.rvdw),
218                        "With Coulomb and LJ PME, rcoulomb should be equal to rvdw");
219
220     pme_lb = new pme_load_balancing_t;
221
222     pme_lb->bSepPMERanks = !thisRankHasDuty(cr, DUTY_PME);
223
224     /* Initially we turn on balancing directly on based on PP/PME imbalance */
225     pme_lb->bTriggerOnDLB = FALSE;
226
227     /* Any number of stages >= 2 is supported */
228     pme_lb->nstage = 2;
229
230     pme_lb->cutoff_scheme = ir.cutoff_scheme;
231
232     pme_lb->rbufOuter_coulomb = nbv.pairlistOuterRadius() - ic.rcoulomb;
233     pme_lb->rbufOuter_vdw     = nbv.pairlistOuterRadius() - ic.rvdw;
234     pme_lb->rbufInner_coulomb = nbv.pairlistInnerRadius() - ic.rcoulomb;
235     pme_lb->rbufInner_vdw     = nbv.pairlistInnerRadius() - ic.rvdw;
236
237     /* Scale box with Ewald wall factor; note that we pmedata->boxScaler
238      * can't always usedd as it's not available with separate PME ranks.
239      */
240     EwaldBoxZScaler boxScaler(ir);
241     boxScaler.scaleBox(box, pme_lb->box_start);
242
243     pme_lb->setup.resize(1);
244
245     pme_lb->rcut_vdw           = ic.rvdw;
246     pme_lb->rcut_coulomb_start = ir.rcoulomb;
247
248     pme_lb->cur                    = 0;
249     pme_lb->setup[0].rcut_coulomb  = ic.rcoulomb;
250     pme_lb->setup[0].rlistOuter    = nbv.pairlistOuterRadius();
251     pme_lb->setup[0].rlistInner    = nbv.pairlistInnerRadius();
252     pme_lb->setup[0].grid[XX]      = ir.nkx;
253     pme_lb->setup[0].grid[YY]      = ir.nky;
254     pme_lb->setup[0].grid[ZZ]      = ir.nkz;
255     pme_lb->setup[0].ewaldcoeff_q  = ic.ewaldcoeff_q;
256     pme_lb->setup[0].ewaldcoeff_lj = ic.ewaldcoeff_lj;
257
258     if (!pme_lb->bSepPMERanks)
259     {
260         GMX_RELEASE_ASSERT(pmedata,
261                            "On ranks doing both PP and PME we need a valid pmedata object");
262         pme_lb->setup[0].pmedata = pmedata;
263     }
264
265     spm = 0;
266     for (d = 0; d < DIM; d++)
267     {
268         sp = norm(pme_lb->box_start[d]) / pme_lb->setup[0].grid[d];
269         if (sp > spm)
270         {
271             spm = sp;
272         }
273     }
274     pme_lb->setup[0].spacing = spm;
275
276     if (ir.fourier_spacing > 0)
277     {
278         pme_lb->cut_spacing = ir.rcoulomb / ir.fourier_spacing;
279     }
280     else
281     {
282         pme_lb->cut_spacing = ir.rcoulomb / pme_lb->setup[0].spacing;
283     }
284
285     pme_lb->stage = 0;
286
287     pme_lb->fastest     = 0;
288     pme_lb->lower_limit = 0;
289     pme_lb->start       = 0;
290     pme_lb->end         = 0;
291     pme_lb->elimited    = epmelblimNO;
292
293     pme_lb->cycles_n  = 0;
294     pme_lb->cycles_c  = 0;
295     pme_lb->startTime = gmx_gettime();
296
297     if (!wallcycle_have_counter())
298     {
299         GMX_LOG(mdlog.warning)
300                 .asParagraph()
301                 .appendText(
302                         "NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use "
303                         "PME-PP balancing.");
304     }
305
306     /* Tune with GPUs and/or separate PME ranks.
307      * When running only on a CPU without PME ranks, PME tuning will only help
308      * with small numbers of atoms in the cut-off sphere.
309      */
310     pme_lb->bActive = (wallcycle_have_counter() && (bUseGPU || pme_lb->bSepPMERanks));
311
312     /* With GPUs and no separate PME ranks we can't measure the PP/PME
313      * imbalance, so we start balancing right away.
314      * Otherwise we only start balancing after we observe imbalance.
315      */
316     pme_lb->bBalance = (pme_lb->bActive && (bUseGPU && !pme_lb->bSepPMERanks));
317
318     pme_lb->step_rel_stop = PMETunePeriod * ir.nstlist;
319
320     /* Delay DD load balancing when GPUs are used */
321     if (pme_lb->bActive && DOMAINDECOMP(cr) && cr->dd->nnodes > 1 && bUseGPU)
322     {
323         /* Lock DLB=auto to off (does nothing when DLB=yes/no.
324          * With GPUs and separate PME nodes, we want to first
325          * do PME tuning without DLB, since DLB might limit
326          * the cut-off, which never improves performance.
327          * We allow for DLB + PME tuning after a first round of tuning.
328          */
329         dd_dlb_lock(cr->dd);
330         if (dd_dlb_is_locked(cr->dd))
331         {
332             GMX_LOG(mdlog.warning)
333                     .asParagraph()
334                     .appendText("NOTE: DLB will not turn on during the first phase of PME tuning");
335         }
336     }
337
338     *pme_lb_p = pme_lb;
339 }
340
341 /*! \brief Try to increase the cutoff during load balancing */
342 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t* pme_lb, int pme_order, const gmx_domdec_t* dd)
343 {
344     real fac, sp;
345     real tmpr_coulomb, tmpr_vdw;
346     int  d;
347     bool grid_ok;
348
349     /* Try to add a new setup with next larger cut-off to the list */
350     pme_setup_t set;
351
352     set.pmedata = nullptr;
353
354     NumPmeDomains numPmeDomains = getNumPmeDomains(dd);
355
356     fac = 1;
357     do
358     {
359         /* Avoid infinite while loop, which can occur at the minimum grid size.
360          * Note that in practice load balancing will stop before this point.
361          * The factor 2.1 allows for the extreme case in which only grids
362          * of powers of 2 are allowed (the current code supports more grids).
363          */
364         if (fac > 2.1)
365         {
366             return FALSE;
367         }
368
369         fac *= 1.01;
370         clear_ivec(set.grid);
371         sp = calcFftGrid(nullptr, pme_lb->box_start, fac * pme_lb->setup[pme_lb->cur].spacing,
372                          minimalPmeGridSize(pme_order), &set.grid[XX], &set.grid[YY], &set.grid[ZZ]);
373
374         /* As here we can't easily check if one of the PME ranks
375          * uses threading, we do a conservative grid check.
376          * This means we can't use pme_order or less grid lines
377          * per PME rank along x, which is not a strong restriction.
378          */
379         grid_ok = gmx_pme_check_restrictions(pme_order, set.grid[XX], set.grid[YY], set.grid[ZZ],
380                                              numPmeDomains.x, true, false);
381     } while (sp <= 1.001 * pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
382
383     set.rcut_coulomb = pme_lb->cut_spacing * sp;
384     if (set.rcut_coulomb < pme_lb->rcut_coulomb_start)
385     {
386         /* This is unlikely, but can happen when e.g. continuing from
387          * a checkpoint after equilibration where the box shrank a lot.
388          * We want to avoid rcoulomb getting smaller than rvdw
389          * and there might be more issues with decreasing rcoulomb.
390          */
391         set.rcut_coulomb = pme_lb->rcut_coulomb_start;
392     }
393
394     if (pme_lb->cutoff_scheme == ecutsVERLET)
395     {
396         /* Never decrease the Coulomb and VdW list buffers */
397         set.rlistOuter = std::max(set.rcut_coulomb + pme_lb->rbufOuter_coulomb,
398                                   pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw);
399         set.rlistInner = std::max(set.rcut_coulomb + pme_lb->rbufInner_coulomb,
400                                   pme_lb->rcut_vdw + pme_lb->rbufInner_vdw);
401     }
402     else
403     {
404         /* TODO Remove these lines and pme_lb->cutoff_scheme */
405         tmpr_coulomb = set.rcut_coulomb + pme_lb->rbufOuter_coulomb;
406         tmpr_vdw     = pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw;
407         /* Two (known) bugs with cutoff-scheme=group here:
408          * - This modification of rlist results in incorrect DD comunication.
409          * - We should set fr->bTwinRange = (fr->rlistlong > fr->rlist).
410          */
411         set.rlistOuter = std::min(tmpr_coulomb, tmpr_vdw);
412         set.rlistInner = set.rlistOuter;
413     }
414
415     set.spacing = sp;
416     /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
417     set.grid_efficiency = 1;
418     for (d = 0; d < DIM; d++)
419     {
420         set.grid_efficiency *= (set.grid[d] * sp) / norm(pme_lb->box_start[d]);
421     }
422     /* The Ewald coefficient is inversly proportional to the cut-off */
423     set.ewaldcoeff_q = pme_lb->setup[0].ewaldcoeff_q * pme_lb->setup[0].rcut_coulomb / set.rcut_coulomb;
424     /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
425     set.ewaldcoeff_lj = pme_lb->setup[0].ewaldcoeff_lj * pme_lb->setup[0].rcut_coulomb / set.rcut_coulomb;
426
427     set.count  = 0;
428     set.cycles = 0;
429
430     if (debug)
431     {
432         fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n", set.grid[XX],
433                 set.grid[YY], set.grid[ZZ], set.rcut_coulomb);
434     }
435     pme_lb->setup.push_back(set);
436     return TRUE;
437 }
438
439 /*! \brief Print the PME grid */
440 static void print_grid(FILE* fp_err, FILE* fp_log, const char* pre, const char* desc, const pme_setup_t* set, double cycles)
441 {
442     auto buf = gmx::formatString("%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f", pre, desc,
443                                  set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
444     if (cycles >= 0)
445     {
446         buf += gmx::formatString(": %.1f M-cycles", cycles * 1e-6);
447     }
448     if (fp_err != nullptr)
449     {
450         fprintf(fp_err, "\r%s\n", buf.c_str());
451         fflush(fp_err);
452     }
453     if (fp_log != nullptr)
454     {
455         fprintf(fp_log, "%s\n", buf.c_str());
456     }
457 }
458
459 /*! \brief Return the index of the last setup used in PME load balancing */
460 static int pme_loadbal_end(pme_load_balancing_t* pme_lb)
461 {
462     /* In the initial stage only n is set; end is not set yet */
463     if (pme_lb->end > 0)
464     {
465         return pme_lb->end;
466     }
467     else
468     {
469         return pme_lb->setup.size();
470     }
471 }
472
473 /*! \brief Print descriptive string about what limits PME load balancing */
474 static void print_loadbal_limited(FILE* fp_err, FILE* fp_log, int64_t step, pme_load_balancing_t* pme_lb)
475 {
476     auto buf = gmx::formatString(
477             "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
478             gmx::int64ToString(step).c_str(), pmelblim_str[pme_lb->elimited],
479             pme_lb->setup[pme_loadbal_end(pme_lb) - 1].rcut_coulomb);
480     if (fp_err != nullptr)
481     {
482         fprintf(fp_err, "\r%s\n", buf.c_str());
483         fflush(fp_err);
484     }
485     if (fp_log != nullptr)
486     {
487         fprintf(fp_log, "%s\n", buf.c_str());
488     }
489 }
490
491 /*! \brief Switch load balancing to stage 1
492  *
493  * In this stage, only reasonably fast setups are run again. */
494 static void switch_to_stage1(pme_load_balancing_t* pme_lb)
495 {
496     /* Increase start until we find a setup that is not slower than
497      * maxRelativeSlowdownAccepted times the fastest setup.
498      */
499     pme_lb->start = pme_lb->lower_limit;
500     while (pme_lb->start + 1 < gmx::ssize(pme_lb->setup)
501            && (pme_lb->setup[pme_lb->start].count == 0
502                || pme_lb->setup[pme_lb->start].cycles
503                           > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted))
504     {
505         pme_lb->start++;
506     }
507     /* While increasing start, we might have skipped setups that we did not
508      * time during stage 0. We want to extend the range for stage 1 to include
509      * any skipped setups that lie between setups that were measured to be
510      * acceptably fast and too slow.
511      */
512     while (pme_lb->start > pme_lb->lower_limit && pme_lb->setup[pme_lb->start - 1].count == 0)
513     {
514         pme_lb->start--;
515     }
516
517     /* Decrease end only with setups that we timed and that are slow. */
518     pme_lb->end = pme_lb->setup.size();
519     if (pme_lb->setup[pme_lb->end - 1].count > 0
520         && pme_lb->setup[pme_lb->end - 1].cycles
521                    > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted)
522     {
523         pme_lb->end--;
524     }
525
526     pme_lb->stage = 1;
527
528     /* Next we want to choose setup pme_lb->end-1, but as we will decrease
529      * pme_lb->cur by one right after returning, we set cur to end.
530      */
531     pme_lb->cur = pme_lb->end;
532 }
533
534 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
535  *
536  * The adjustment is done to generate a different non-bonded PP and PME load.
537  * With separate PME ranks (PP and PME on different processes) or with
538  * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
539  * and changing the load will affect the load balance and performance.
540  * The total time for a set of integration steps is monitored and a range
541  * of grid/cut-off setups is scanned. After calling pme_load_balance many
542  * times and acquiring enough statistics, the best performing setup is chosen.
543  * Here we try to take into account fluctuations and changes due to external
544  * factors as well as DD load balancing.
545  */
546 static void pme_load_balance(pme_load_balancing_t*          pme_lb,
547                              t_commrec*                     cr,
548                              FILE*                          fp_err,
549                              FILE*                          fp_log,
550                              const gmx::MDLogger&           mdlog,
551                              const t_inputrec&              ir,
552                              const matrix                   box,
553                              gmx::ArrayRef<const gmx::RVec> x,
554                              double                         cycles,
555                              interaction_const_t*           ic,
556                              struct nonbonded_verlet_t*     nbv,
557                              struct gmx_pme_t**             pmedata,
558                              int64_t                        step)
559 {
560     gmx_bool     OK;
561     pme_setup_t* set;
562     double       cycles_fast;
563     char         buf[STRLEN], sbuf[22];
564
565     if (PAR(cr))
566     {
567         gmx_sumd(1, &cycles, cr);
568         cycles /= cr->nnodes;
569     }
570
571     set = &pme_lb->setup[pme_lb->cur];
572     set->count++;
573
574     /* Skip the first c_numPostSwitchTuningIntervalSkip cycles because the first step
575      * after a switch is much slower due to allocation and/or caching effects.
576      */
577     if (set->count % (c_numPostSwitchTuningIntervalSkip + 1) != 0)
578     {
579         return;
580     }
581
582     sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
583     print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
584
585     GMX_RELEASE_ASSERT(set->count > c_numPostSwitchTuningIntervalSkip, "We should skip cycles");
586     if (set->count == (c_numPostSwitchTuningIntervalSkip + 1))
587     {
588         set->cycles = cycles;
589     }
590     else
591     {
592         if (cycles * maxFluctuationAccepted < set->cycles && pme_lb->stage == pme_lb->nstage - 1)
593         {
594             /* The performance went up a lot (due to e.g. DD load balancing).
595              * Add a stage, keep the minima, but rescan all setups.
596              */
597             pme_lb->nstage++;
598
599             if (debug)
600             {
601                 fprintf(debug,
602                         "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this "
603                         "is more than %f\n"
604                         "Increased the number stages to %d"
605                         " and ignoring the previous performance\n",
606                         set->grid[XX], set->grid[YY], set->grid[ZZ], set->cycles * 1e-6,
607                         cycles * 1e-6, maxFluctuationAccepted, pme_lb->nstage);
608             }
609         }
610         set->cycles = std::min(set->cycles, cycles);
611     }
612
613     if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
614     {
615         pme_lb->fastest = pme_lb->cur;
616
617         if (DOMAINDECOMP(cr))
618         {
619             /* We found a new fastest setting, ensure that with subsequent
620              * shorter cut-off's the dynamic load balancing does not make
621              * the use of the current cut-off impossible. This solution is
622              * a trade-off, as the PME load balancing and DD domain size
623              * load balancing can interact in complex ways.
624              * With the Verlet kernels, DD load imbalance will usually be
625              * mainly due to bonded interaction imbalance, which will often
626              * quickly push the domain boundaries beyond the limit for the
627              * optimal, PME load balanced, cut-off. But it could be that
628              * better overal performance can be obtained with a slightly
629              * shorter cut-off and better DD load balancing.
630              */
631             set_dd_dlb_max_cutoff(cr, pme_lb->setup[pme_lb->fastest].rlistOuter);
632         }
633     }
634     cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
635
636     /* Check in stage 0 if we should stop scanning grids.
637      * Stop when the time is more than maxRelativeSlowDownAccepted longer than the fastest.
638      */
639     if (pme_lb->stage == 0 && pme_lb->cur > 0
640         && cycles > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted)
641     {
642         pme_lb->setup.resize(pme_lb->cur + 1);
643         /* Done with scanning, go to stage 1 */
644         switch_to_stage1(pme_lb);
645     }
646
647     if (pme_lb->stage == 0)
648     {
649         int gridsize_start;
650
651         gridsize_start = set->grid[XX] * set->grid[YY] * set->grid[ZZ];
652
653         do
654         {
655             if (pme_lb->cur + 1 < gmx::ssize(pme_lb->setup))
656             {
657                 /* We had already generated the next setup */
658                 OK = TRUE;
659             }
660             else
661             {
662                 /* Find the next setup */
663                 OK = pme_loadbal_increase_cutoff(pme_lb, ir.pme_order, cr->dd);
664
665                 if (!OK)
666                 {
667                     pme_lb->elimited = epmelblimPMEGRID;
668                 }
669             }
670
671             if (OK
672                 && pme_lb->setup[pme_lb->cur + 1].spacing > c_maxSpacingScaling * pme_lb->setup[0].spacing)
673             {
674                 OK               = FALSE;
675                 pme_lb->elimited = epmelblimMAXSCALING;
676             }
677
678             if (OK && ir.ePBC != epbcNONE)
679             {
680                 OK = (gmx::square(pme_lb->setup[pme_lb->cur + 1].rlistOuter) <= max_cutoff2(ir.ePBC, box));
681                 if (!OK)
682                 {
683                     pme_lb->elimited = epmelblimBOX;
684                 }
685             }
686
687             if (OK)
688             {
689                 pme_lb->cur++;
690
691                 if (DOMAINDECOMP(cr))
692                 {
693                     OK = change_dd_cutoff(cr, box, x, pme_lb->setup[pme_lb->cur].rlistOuter);
694                     if (!OK)
695                     {
696                         /* Failed: do not use this setup */
697                         pme_lb->cur--;
698                         pme_lb->elimited = epmelblimDD;
699                     }
700                 }
701             }
702             if (!OK)
703             {
704                 /* We hit the upper limit for the cut-off,
705                  * the setup should not go further than cur.
706                  */
707                 pme_lb->setup.resize(pme_lb->cur + 1);
708                 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
709                 /* Switch to the next stage */
710                 switch_to_stage1(pme_lb);
711             }
712         } while (OK
713                  && !(pme_lb->setup[pme_lb->cur].grid[XX] * pme_lb->setup[pme_lb->cur].grid[YY]
714                                       * pme_lb->setup[pme_lb->cur].grid[ZZ]
715                               < gridsize_start * gridpointsScaleFactor
716                       && pme_lb->setup[pme_lb->cur].grid_efficiency
717                                  < pme_lb->setup[pme_lb->cur - 1].grid_efficiency * relativeEfficiencyFactor));
718     }
719
720     if (pme_lb->stage > 0 && pme_lb->end == 1)
721     {
722         pme_lb->cur   = pme_lb->lower_limit;
723         pme_lb->stage = pme_lb->nstage;
724     }
725     else if (pme_lb->stage > 0 && pme_lb->end > 1)
726     {
727         /* If stage = nstage-1:
728          *   scan over all setups, rerunning only those setups
729          *   which are not much slower than the fastest
730          * else:
731          *   use the next setup
732          * Note that we loop backward to minimize the risk of the cut-off
733          * getting limited by DD DLB, since the DLB cut-off limit is set
734          * to the fastest PME setup.
735          */
736         do
737         {
738             if (pme_lb->cur > pme_lb->start)
739             {
740                 pme_lb->cur--;
741             }
742             else
743             {
744                 pme_lb->stage++;
745
746                 pme_lb->cur = pme_lb->end - 1;
747             }
748         } while (pme_lb->stage == pme_lb->nstage - 1 && pme_lb->setup[pme_lb->cur].count > 0
749                  && pme_lb->setup[pme_lb->cur].cycles > cycles_fast * maxRelativeSlowdownAccepted);
750
751         if (pme_lb->stage == pme_lb->nstage)
752         {
753             /* We are done optimizing, use the fastest setup we found */
754             pme_lb->cur = pme_lb->fastest;
755         }
756     }
757
758     if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
759     {
760         OK = change_dd_cutoff(cr, box, x, pme_lb->setup[pme_lb->cur].rlistOuter);
761         if (!OK)
762         {
763             /* For some reason the chosen cut-off is incompatible with DD.
764              * We should continue scanning a more limited range of cut-off's.
765              */
766             if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
767             {
768                 /* stage=nstage says we're finished, but we should continue
769                  * balancing, so we set back stage which was just incremented.
770                  */
771                 pme_lb->stage--;
772             }
773             if (pme_lb->cur <= pme_lb->fastest)
774             {
775                 /* This should not happen, as we set limits on the DLB bounds.
776                  * But we implement a complete failsafe solution anyhow.
777                  */
778                 GMX_LOG(mdlog.warning)
779                         .asParagraph()
780                         .appendTextFormatted(
781                                 "The fastest PP/PME load balancing setting (cutoff %.3d nm) is no "
782                                 "longer available due to DD DLB or box size limitations",
783                                 pme_lb->fastest);
784                 pme_lb->fastest = pme_lb->lower_limit;
785                 pme_lb->start   = pme_lb->lower_limit;
786             }
787             /* Limit the range to below the current cut-off, scan from start */
788             pme_lb->end      = pme_lb->cur;
789             pme_lb->cur      = pme_lb->start;
790             pme_lb->elimited = epmelblimDD;
791             print_loadbal_limited(fp_err, fp_log, step, pme_lb);
792         }
793     }
794
795     /* Change the Coulomb cut-off and the PME grid */
796
797     set = &pme_lb->setup[pme_lb->cur];
798
799     ic->rcoulomb = set->rcut_coulomb;
800     nbv->changePairlistRadii(set->rlistOuter, set->rlistInner);
801     ic->ewaldcoeff_q = set->ewaldcoeff_q;
802     /* TODO: centralize the code that sets the potentials shifts */
803     if (ic->coulomb_modifier == eintmodPOTSHIFT)
804     {
805         GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
806         ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
807     }
808     if (EVDW_PME(ic->vdwtype))
809     {
810         /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
811         ic->rvdw          = set->rcut_coulomb;
812         ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
813         if (ic->vdw_modifier == eintmodPOTSHIFT)
814         {
815             real crc2;
816
817             ic->dispersion_shift.cpot = -1.0 / gmx::power6(static_cast<double>(ic->rvdw));
818             ic->repulsion_shift.cpot  = -1.0 / gmx::power12(static_cast<double>(ic->rvdw));
819             crc2                      = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
820             ic->sh_lj_ewald =
821                     (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
822         }
823     }
824
825     /* We always re-initialize the tables whether they are used or not */
826     init_interaction_const_tables(nullptr, ic);
827
828     Nbnxm::gpu_pme_loadbal_update_param(nbv, ic);
829
830     if (!pme_lb->bSepPMERanks)
831     {
832         /* FIXME:
833          * CPU PME keeps a list of allocated pmedata's, that's why pme_lb->setup[pme_lb->cur].pmedata is not always nullptr.
834          * GPU PME, however, currently needs the gmx_pme_reinit always called on load balancing
835          * (pme_gpu_reinit might be not sufficiently decoupled from gmx_pme_init).
836          * This can lead to a lot of reallocations for PME GPU.
837          * Would be nicer if the allocated grid list was hidden within a single pmedata structure.
838          */
839         if ((pme_lb->setup[pme_lb->cur].pmedata == nullptr)
840             || pme_gpu_task_enabled(pme_lb->setup[pme_lb->cur].pmedata))
841         {
842             /* Generate a new PME data structure,
843              * copying part of the old pointers.
844              */
845             gmx_pme_reinit(&set->pmedata, cr, pme_lb->setup[0].pmedata, &ir, set->grid,
846                            set->ewaldcoeff_q, set->ewaldcoeff_lj);
847         }
848         *pmedata = set->pmedata;
849     }
850     else
851     {
852         /* Tell our PME-only rank to switch grid */
853         gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
854     }
855
856     if (debug)
857     {
858         print_grid(nullptr, debug, "", "switched to", set, -1);
859     }
860
861     if (pme_lb->stage == pme_lb->nstage)
862     {
863         print_grid(fp_err, fp_log, "", "optimal", set, -1);
864     }
865 }
866
867 /*! \brief Prepare for another round of PME load balancing
868  *
869  * \param[in,out] pme_lb       Pointer to PME load balancing struct
870  * \param[in]     bDlbUnlocked TRUE is DLB was locked and is now unlocked
871  *
872  * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
873  * the PP/PME balance might change and re-balancing can improve performance.
874  * This function adds 2 stages and adjusts the considered setup range.
875  */
876 static void continue_pme_loadbal(pme_load_balancing_t* pme_lb, gmx_bool bDlbUnlocked)
877 {
878     /* Add 2 tuning stages, keep the detected end of the setup range */
879     pme_lb->nstage += 2;
880     if (bDlbUnlocked && pme_lb->bSepPMERanks)
881     {
882         /* With separate PME ranks, DLB should always lower the PP load and
883          * can only increase the PME load (more communication and imbalance),
884          * so we only need to scan longer cut-off's.
885          */
886         pme_lb->lower_limit = pme_lb->cur;
887     }
888     pme_lb->start = pme_lb->lower_limit;
889 }
890
891 void pme_loadbal_do(pme_load_balancing_t*          pme_lb,
892                     t_commrec*                     cr,
893                     FILE*                          fp_err,
894                     FILE*                          fp_log,
895                     const gmx::MDLogger&           mdlog,
896                     const t_inputrec&              ir,
897                     t_forcerec*                    fr,
898                     const matrix                   box,
899                     gmx::ArrayRef<const gmx::RVec> x,
900                     gmx_wallcycle_t                wcycle,
901                     int64_t                        step,
902                     int64_t                        step_rel,
903                     gmx_bool*                      bPrinting,
904                     bool                           useGpuPmePpCommunication)
905 {
906     int    n_prev;
907     double cycles_prev;
908
909     assert(pme_lb != nullptr);
910
911     if (!pme_lb->bActive)
912     {
913         return;
914     }
915
916     n_prev      = pme_lb->cycles_n;
917     cycles_prev = pme_lb->cycles_c;
918     wallcycle_get(wcycle, ewcSTEP, &pme_lb->cycles_n, &pme_lb->cycles_c);
919
920     /* Before the first step we haven't done any steps yet.
921      * Also handle cases where ir.init_step % ir.nstlist != 0.
922      * We also want to skip a number of steps and seconds while
923      * the CPU and GPU, when used, performance stabilizes.
924      */
925     if (pme_lb->cycles_n == 0 || step_rel < c_numFirstTuningIntervalSkip * ir.nstlist
926         || gmx_gettime() - pme_lb->startTime < c_startupTimeDelay)
927     {
928         *bPrinting = FALSE;
929
930         return;
931     }
932     /* Sanity check, we expect nstlist cycle counts */
933     if (pme_lb->cycles_n - n_prev != ir.nstlist)
934     {
935         /* We could return here, but it's safer to issue an error and quit */
936         gmx_incons("pme_loadbal_do called at an interval != nstlist");
937     }
938
939     /* PME grid + cut-off optimization with GPUs or PME ranks */
940     if (!pme_lb->bBalance && pme_lb->bSepPMERanks)
941     {
942         if (pme_lb->bTriggerOnDLB)
943         {
944             pme_lb->bBalance = dd_dlb_is_on(cr->dd);
945         }
946         /* We should ignore the first timing to avoid timing allocation
947          * overhead. And since the PME load balancing is called just
948          * before DD repartitioning, the ratio returned by dd_pme_f_ratio
949          * is not over the last nstlist steps, but the nstlist steps before
950          * that. So the first useful ratio is available at step_rel=3*nstlist.
951          */
952         else if (step_rel >= c_numFirstTuningIntervalSkipWithSepPme * ir.nstlist)
953         {
954             if (DDMASTER(cr->dd))
955             {
956                 /* If PME rank load is too high, start tuning. If
957                    PME-PP direct GPU communication is active,
958                    unconditionally start tuning since ratio will be
959                    unreliable due to CPU-GPU asynchronicity in codepath */
960                 pme_lb->bBalance = useGpuPmePpCommunication
961                                            ? true
962                                            : (dd_pme_f_ratio(cr->dd) >= loadBalanceTriggerFactor);
963             }
964             dd_bcast(cr->dd, sizeof(gmx_bool), &pme_lb->bBalance);
965         }
966
967         pme_lb->bActive = (pme_lb->bBalance || step_rel <= pme_lb->step_rel_stop);
968     }
969
970     /* The location in the code of this balancing termination is strange.
971      * You would expect to have it after the call to pme_load_balance()
972      * below, since there pme_lb->stage is updated.
973      * But when terminating directly after deciding on and selecting the
974      * optimal setup, DLB will turn on right away if it was locked before.
975      * This might be due to PME reinitialization. So we check stage here
976      * to allow for another nstlist steps with DLB locked to stabilize
977      * the performance.
978      */
979     if (pme_lb->bBalance && pme_lb->stage == pme_lb->nstage)
980     {
981         pme_lb->bBalance = FALSE;
982
983         if (DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
984         {
985             /* Unlock the DLB=auto, DLB is allowed to activate */
986             dd_dlb_unlock(cr->dd);
987             GMX_LOG(mdlog.warning)
988                     .asParagraph()
989                     .appendText("NOTE: DLB can now turn on, when beneficial");
990
991             /* We don't deactivate the tuning yet, since we will balance again
992              * after DLB gets turned on, if it does within PMETune_period.
993              */
994             continue_pme_loadbal(pme_lb, TRUE);
995             pme_lb->bTriggerOnDLB = TRUE;
996             pme_lb->step_rel_stop = step_rel + PMETunePeriod * ir.nstlist;
997         }
998         else
999         {
1000             /* We're completely done with PME tuning */
1001             pme_lb->bActive = FALSE;
1002         }
1003
1004         if (DOMAINDECOMP(cr))
1005         {
1006             /* Set the cut-off limit to the final selected cut-off,
1007              * so we don't have artificial DLB limits.
1008              * This also ensures that we won't disable the currently
1009              * optimal setting during a second round of PME balancing.
1010              */
1011             set_dd_dlb_max_cutoff(cr, fr->nbv->pairlistOuterRadius());
1012         }
1013     }
1014
1015     if (pme_lb->bBalance)
1016     {
1017         /* We might not have collected nstlist steps in cycles yet,
1018          * since init_step might not be a multiple of nstlist,
1019          * but the first data collected is skipped anyhow.
1020          */
1021         pme_load_balance(pme_lb, cr, fp_err, fp_log, mdlog, ir, box, x,
1022                          pme_lb->cycles_c - cycles_prev, fr->ic, fr->nbv.get(), &fr->pmedata, step);
1023
1024         /* Update deprecated rlist in forcerec to stay in sync with fr->nbv */
1025         fr->rlist = fr->nbv->pairlistOuterRadius();
1026
1027         if (ir.eDispCorr != edispcNO)
1028         {
1029             fr->dispersionCorrection->setParameters(*fr->ic);
1030         }
1031     }
1032
1033     if (!pme_lb->bBalance && (!pme_lb->bSepPMERanks || step_rel > pme_lb->step_rel_stop))
1034     {
1035         /* We have just deactivated the balancing and we're not measuring PP/PME
1036          * imbalance during the first steps of the run: deactivate the tuning.
1037          */
1038         pme_lb->bActive = FALSE;
1039     }
1040
1041     if (!(pme_lb->bActive) && DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
1042     {
1043         /* Make sure DLB is allowed when we deactivate PME tuning */
1044         dd_dlb_unlock(cr->dd);
1045         GMX_LOG(mdlog.warning)
1046                 .asParagraph()
1047                 .appendText("NOTE: DLB can now turn on, when beneficial");
1048     }
1049
1050     *bPrinting = pme_lb->bBalance;
1051 }
1052
1053 /*! \brief Return product of the number of PME grid points in each dimension */
1054 static int pme_grid_points(const pme_setup_t* setup)
1055 {
1056     return setup->grid[XX] * setup->grid[YY] * setup->grid[ZZ];
1057 }
1058
1059 /*! \brief Print one load-balancing setting */
1060 static void print_pme_loadbal_setting(FILE* fplog, const char* name, const pme_setup_t* setup)
1061 {
1062     fprintf(fplog, "   %-7s %6.3f nm %6.3f nm     %3d %3d %3d   %5.3f nm  %5.3f nm\n", name,
1063             setup->rcut_coulomb, setup->rlistInner, setup->grid[XX], setup->grid[YY],
1064             setup->grid[ZZ], setup->spacing, 1 / setup->ewaldcoeff_q);
1065 }
1066
1067 /*! \brief Print all load-balancing settings */
1068 static void print_pme_loadbal_settings(pme_load_balancing_t* pme_lb,
1069                                        FILE*                 fplog,
1070                                        const gmx::MDLogger&  mdlog,
1071                                        gmx_bool              bNonBondedOnGPU)
1072 {
1073     double pp_ratio, grid_ratio;
1074     real   pp_ratio_temporary;
1075
1076     pp_ratio_temporary = pme_lb->setup[pme_lb->cur].rlistInner / pme_lb->setup[0].rlistInner;
1077     pp_ratio           = gmx::power3(pp_ratio_temporary);
1078     grid_ratio         = pme_grid_points(&pme_lb->setup[pme_lb->cur])
1079                  / static_cast<double>(pme_grid_points(&pme_lb->setup[0]));
1080
1081     fprintf(fplog, "\n");
1082     fprintf(fplog, "       P P   -   P M E   L O A D   B A L A N C I N G\n");
1083     fprintf(fplog, "\n");
1084     /* Here we only warn when the optimal setting is the last one */
1085     if (pme_lb->elimited != epmelblimNO && pme_lb->cur == pme_loadbal_end(pme_lb) - 1)
1086     {
1087         fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
1088                 pmelblim_str[pme_lb->elimited]);
1089         fprintf(fplog, "       you might not have reached a good load balance.\n");
1090         if (pme_lb->elimited == epmelblimDD)
1091         {
1092             fprintf(fplog, "       Try different mdrun -dd settings or lower the -dds value.\n");
1093         }
1094         fprintf(fplog, "\n");
1095     }
1096     fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
1097     fprintf(fplog, "           particle-particle                    PME\n");
1098     fprintf(fplog, "            rcoulomb  rlist            grid      spacing   1/beta\n");
1099     print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
1100     print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
1101     fprintf(fplog, " cost-ratio           %4.2f             %4.2f\n", pp_ratio, grid_ratio);
1102     fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
1103
1104     if (pp_ratio > 1.5 && !bNonBondedOnGPU)
1105     {
1106         GMX_LOG(mdlog.warning)
1107                 .asParagraph()
1108                 .appendText(
1109                         "NOTE: PME load balancing increased the non-bonded workload by more than "
1110                         "50%.\n"
1111                         "      For better performance, use (more) PME ranks (mdrun -npme),\n"
1112                         "      or if you are beyond the scaling limit, use fewer total ranks (or "
1113                         "nodes).");
1114     }
1115     else
1116     {
1117         fprintf(fplog, "\n");
1118     }
1119 }
1120
1121 void pme_loadbal_done(pme_load_balancing_t* pme_lb, FILE* fplog, const gmx::MDLogger& mdlog, gmx_bool bNonBondedOnGPU)
1122 {
1123     if (fplog != nullptr && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
1124     {
1125         print_pme_loadbal_settings(pme_lb, fplog, mdlog, bNonBondedOnGPU);
1126     }
1127
1128     delete pme_lb;
1129 }