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