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