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