2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief This file contains function definitions necessary for
39 * managing automatic load balance of PME calculations (Coulomb and
42 * \author Berk Hess <hess@kth.se>
43 * \ingroup module_ewald
47 #include "gromacs/utility/enumerationhelpers.h"
48 #include "pme_load_balancing.h"
55 #include "gromacs/domdec/dlb.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/partition.h"
60 #include "gromacs/ewald/ewald_utils.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fft/calcgrid.h"
63 #include "gromacs/gmxlib/network.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/vec.h"
66 #include "gromacs/mdlib/dispersioncorrection.h"
67 #include "gromacs/mdlib/forcerec.h"
68 #include "gromacs/mdtypes/commrec.h"
69 #include "gromacs/mdtypes/forcerec.h"
70 #include "gromacs/mdtypes/inputrec.h"
71 #include "gromacs/mdtypes/interaction_const.h"
72 #include "gromacs/mdtypes/md_enums.h"
73 #include "gromacs/mdtypes/state.h"
74 #include "gromacs/nbnxm/gpu_data_mgmt.h"
75 #include "gromacs/nbnxm/nbnxm.h"
76 #include "gromacs/pbcutil/pbc.h"
77 #include "gromacs/timing/wallcycle.h"
78 #include "gromacs/timing/walltime_accounting.h"
79 #include "gromacs/utility/cstringutil.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/gmxassert.h"
82 #include "gromacs/utility/logger.h"
83 #include "gromacs/utility/smalloc.h"
84 #include "gromacs/utility/strconvert.h"
86 #include "pme_internal.h"
89 /*! \brief Parameters and settings for one PP-PME setup */
92 real rcut_coulomb; /**< Coulomb cut-off */
93 real rlistOuter; /**< cut-off for the outer pair-list */
94 real rlistInner; /**< cut-off for the inner pair-list */
95 real spacing; /**< (largest) PME grid spacing */
96 ivec grid; /**< the PME grid dimensions */
97 real grid_efficiency; /**< ineffiency factor for non-uniform grids <= 1 */
98 real ewaldcoeff_q; /**< Electrostatic Ewald coefficient */
99 real ewaldcoeff_lj; /**< LJ Ewald coefficient, only for the call to send_switchgrid */
100 struct gmx_pme_t* pmedata; /**< the data structure used in the PME code */
101 int count; /**< number of times this setup has been timed */
102 double cycles; /**< the fastest time for this setup in cycles */
105 /*! \brief After 50 nstlist periods of not observing imbalance: never tune PME */
106 const int PMETunePeriod = 50;
107 /*! \brief Trigger PME load balancing at more than 5% PME overload */
108 const real loadBalanceTriggerFactor = 1.05;
109 /*! \brief Scale the grid by a most at factor 1.7.
111 * This still leaves room for about 4-4.5x decrease in grid spacing while limiting the cases where
112 * large imbalance leads to extreme cutoff scaling for marginal benefits.
114 * This should help to avoid:
115 * - large increase in power consumption for little performance gain
116 * - increasing communication volume
119 const real c_maxSpacingScaling = 1.7;
120 /*! \brief In the initial scan, step by grids that are at least a factor 0.8 coarser */
121 const real gridpointsScaleFactor = 0.8;
122 /*! \brief In the initial scan, try to skip grids with uneven x/y/z spacing,
123 * checking if the "efficiency" is more than 5% worse than the previous grid.
125 const real relativeEfficiencyFactor = 1.05;
126 /*! \brief Rerun until a run is 12% slower setups than the fastest run so far */
127 const real maxRelativeSlowdownAccepted = 1.12;
128 /*! \brief If setups get more than 2% faster, do another round to avoid
129 * choosing a slower setup due to acceleration or fluctuations.
131 const real maxFluctuationAccepted = 1.02;
133 //! \brief Number of nstlist long tuning intervals to skip before starting
134 // load-balancing at the beginning of the run.
135 const int c_numFirstTuningIntervalSkip = 5;
136 //! \brief Number of nstlist long tuning intervals to skip before starting
137 // load-balancing at the beginning of the run with separate PME ranks. */
138 const int c_numFirstTuningIntervalSkipWithSepPme = 3;
139 //! \brief Number of nstlist long tuning intervals to skip after switching to a new setting
141 const int c_numPostSwitchTuningIntervalSkip = 1;
142 //! \brief Number of seconds to delay the tuning at startup to allow processors clocks to ramp up.
143 const double c_startupTimeDelay = 5.0;
145 /*! \brief Enumeration whose values describe the effect limiting the load balancing */
146 enum class PmeLoadBalancingLimit : int
156 /*! \brief Descriptive strings for PmeLoadBalancingLimit \c enumValue */
157 static const char* enumValueToString(PmeLoadBalancingLimit enumValue)
159 constexpr gmx::EnumerationArray<PmeLoadBalancingLimit, const char*> pmeLoadBalancingLimitNames = {
162 "domain decompostion",
163 "PME grid restriction",
164 "maximum allowed grid scaling"
166 return pmeLoadBalancingLimitNames[enumValue];
169 struct pme_load_balancing_t
171 gmx_bool bSepPMERanks; /**< do we have separate PME ranks? */
172 gmx_bool bActive; /**< is PME tuning active? */
173 int64_t step_rel_stop; /**< stop the tuning after this value of step_rel */
174 gmx_bool bTriggerOnDLB; /**< trigger balancing only on DD DLB */
175 gmx_bool bBalance; /**< are we in the balancing phase, i.e. trying different setups? */
176 int nstage; /**< the current maximum number of stages */
177 bool startupTimeDelayElapsed; /**< Has the c_startupTimeDelay elapsed indicating that the balancing can start. */
179 real cut_spacing; /**< the minimum cutoff / PME grid spacing ratio */
180 real rcut_vdw; /**< Vdw cutoff (does not change) */
181 real rcut_coulomb_start; /**< Initial electrostatics cutoff */
182 real rbufOuter_coulomb; /**< the outer pairlist buffer size */
183 real rbufOuter_vdw; /**< the outer pairlist buffer size */
184 real rbufInner_coulomb; /**< the inner pairlist buffer size */
185 real rbufInner_vdw; /**< the inner pairlist buffer size */
186 matrix box_start; /**< the initial simulation box */
187 std::vector<pme_setup_t> setup; /**< the PME+cutoff setups */
188 int cur; /**< the index (in setup) of the current setup */
189 int fastest; /**< index of the fastest setup up till now */
190 int lower_limit; /**< don't go below this setup index */
191 int start; /**< start of setup index range to consider in stage>0 */
192 int end; /**< end of setup index range to consider in stage>0 */
193 PmeLoadBalancingLimit elimited; /**< was the balancing limited, uses enum above */
194 CutoffScheme cutoff_scheme; /**< Verlet or group cut-offs */
196 int stage; /**< the current stage */
198 int cycles_n; /**< step cycle counter cumulative count */
199 double cycles_c; /**< step cycle counter cumulative cycles */
200 double startTime; /**< time stamp when the balancing was started on the master rank (relative to the UNIX epoch start).*/
203 /* TODO The code in this file should call this getter, rather than
204 * read bActive anywhere */
205 bool pme_loadbal_is_active(const pme_load_balancing_t* pme_lb)
207 return pme_lb != nullptr && pme_lb->bActive;
210 // TODO Return a unique_ptr to pme_load_balancing_t
211 void pme_loadbal_init(pme_load_balancing_t** pme_lb_p,
213 const gmx::MDLogger& mdlog,
214 const t_inputrec& ir,
216 const interaction_const_t& ic,
217 const nonbonded_verlet_t& nbv,
222 pme_load_balancing_t* pme_lb;
226 // Note that we don't (yet) support PME load balancing with LJ-PME only.
227 GMX_RELEASE_ASSERT(EEL_PME(ir.coulombtype),
228 "pme_loadbal_init called without PME electrostatics");
229 // To avoid complexity, we require a single cut-off with PME for q+LJ.
230 // This is checked by grompp, but it doesn't hurt to check again.
231 GMX_RELEASE_ASSERT(!(EEL_PME(ir.coulombtype) && EVDW_PME(ir.vdwtype) && ir.rcoulomb != ir.rvdw),
232 "With Coulomb and LJ PME, rcoulomb should be equal to rvdw");
234 pme_lb = new pme_load_balancing_t;
236 pme_lb->bSepPMERanks = !thisRankHasDuty(cr, DUTY_PME);
238 /* Initially we turn on balancing directly on based on PP/PME imbalance */
239 pme_lb->bTriggerOnDLB = FALSE;
241 /* Any number of stages >= 2 is supported */
244 pme_lb->cutoff_scheme = ir.cutoff_scheme;
246 pme_lb->rbufOuter_coulomb = nbv.pairlistOuterRadius() - ic.rcoulomb;
247 pme_lb->rbufOuter_vdw = nbv.pairlistOuterRadius() - ic.rvdw;
248 pme_lb->rbufInner_coulomb = nbv.pairlistInnerRadius() - ic.rcoulomb;
249 pme_lb->rbufInner_vdw = nbv.pairlistInnerRadius() - ic.rvdw;
251 /* Scale box with Ewald wall factor; note that we pmedata->boxScaler
252 * can't always usedd as it's not available with separate PME ranks.
254 EwaldBoxZScaler boxScaler(inputrecPbcXY2Walls(&ir), ir.wall_ewald_zfac);
255 boxScaler.scaleBox(box, pme_lb->box_start);
257 pme_lb->setup.resize(1);
259 pme_lb->rcut_vdw = ic.rvdw;
260 pme_lb->rcut_coulomb_start = ir.rcoulomb;
263 pme_lb->setup[0].rcut_coulomb = ic.rcoulomb;
264 pme_lb->setup[0].rlistOuter = nbv.pairlistOuterRadius();
265 pme_lb->setup[0].rlistInner = nbv.pairlistInnerRadius();
266 pme_lb->setup[0].grid[XX] = ir.nkx;
267 pme_lb->setup[0].grid[YY] = ir.nky;
268 pme_lb->setup[0].grid[ZZ] = ir.nkz;
269 pme_lb->setup[0].ewaldcoeff_q = ic.ewaldcoeff_q;
270 pme_lb->setup[0].ewaldcoeff_lj = ic.ewaldcoeff_lj;
272 if (!pme_lb->bSepPMERanks)
274 GMX_RELEASE_ASSERT(pmedata, "On ranks doing both PP and PME we need a valid pmedata object");
275 pme_lb->setup[0].pmedata = pmedata;
279 for (d = 0; d < DIM; d++)
281 sp = norm(pme_lb->box_start[d]) / pme_lb->setup[0].grid[d];
287 pme_lb->setup[0].spacing = spm;
289 if (ir.fourier_spacing > 0)
291 pme_lb->cut_spacing = ir.rcoulomb / ir.fourier_spacing;
295 pme_lb->cut_spacing = ir.rcoulomb / pme_lb->setup[0].spacing;
301 pme_lb->lower_limit = 0;
304 pme_lb->elimited = PmeLoadBalancingLimit::No;
306 pme_lb->cycles_n = 0;
307 pme_lb->cycles_c = 0;
308 // only master ranks do timing
309 if (!PAR(cr) || (haveDDAtomOrdering(*cr) && DDMASTER(cr->dd)))
311 pme_lb->startTime = gmx_gettime();
314 if (!wallcycle_have_counter())
316 GMX_LOG(mdlog.warning)
319 "NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use "
320 "PME-PP balancing.");
323 /* Tune with GPUs and/or separate PME ranks.
324 * When running only on a CPU without PME ranks, PME tuning will only help
325 * with small numbers of atoms in the cut-off sphere.
327 pme_lb->bActive = (wallcycle_have_counter() && (bUseGPU || pme_lb->bSepPMERanks));
329 /* With GPUs and no separate PME ranks we can't measure the PP/PME
330 * imbalance, so we start balancing right away.
331 * Otherwise we only start balancing after we observe imbalance.
333 pme_lb->bBalance = (pme_lb->bActive && (bUseGPU && !pme_lb->bSepPMERanks));
335 pme_lb->step_rel_stop = PMETunePeriod * ir.nstlist;
337 /* Delay DD load balancing when GPUs are used */
338 if (pme_lb->bActive && haveDDAtomOrdering(*cr) && cr->dd->nnodes > 1 && bUseGPU)
340 /* Lock DLB=auto to off (does nothing when DLB=yes/no.
341 * With GPUs and separate PME nodes, we want to first
342 * do PME tuning without DLB, since DLB might limit
343 * the cut-off, which never improves performance.
344 * We allow for DLB + PME tuning after a first round of tuning.
347 if (dd_dlb_is_locked(cr->dd))
349 GMX_LOG(mdlog.warning)
351 .appendText("NOTE: DLB will not turn on during the first phase of PME tuning");
358 /*! \brief Try to increase the cutoff during load balancing */
359 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t* pme_lb, int pme_order, const gmx_domdec_t* dd)
362 real tmpr_coulomb, tmpr_vdw;
366 /* Try to add a new setup with next larger cut-off to the list */
369 set.pmedata = nullptr;
371 NumPmeDomains numPmeDomains = getNumPmeDomains(dd);
376 /* Avoid infinite while loop, which can occur at the minimum grid size.
377 * Note that in practice load balancing will stop before this point.
378 * The factor 2.1 allows for the extreme case in which only grids
379 * of powers of 2 are allowed (the current code supports more grids).
387 clear_ivec(set.grid);
388 sp = calcFftGrid(nullptr,
390 fac * pme_lb->setup[pme_lb->cur].spacing,
391 minimalPmeGridSize(pme_order),
396 /* As here we can't easily check if one of the PME ranks
397 * uses threading, we do a conservative grid check.
398 * This means we can't use pme_order or less grid lines
399 * per PME rank along x, which is not a strong restriction.
401 grid_ok = gmx_pme_check_restrictions(
402 pme_order, set.grid[XX], set.grid[YY], set.grid[ZZ], numPmeDomains.x, true, false);
403 } while (sp <= 1.001 * pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
405 set.rcut_coulomb = pme_lb->cut_spacing * sp;
406 if (set.rcut_coulomb < pme_lb->rcut_coulomb_start)
408 /* This is unlikely, but can happen when e.g. continuing from
409 * a checkpoint after equilibration where the box shrank a lot.
410 * We want to avoid rcoulomb getting smaller than rvdw
411 * and there might be more issues with decreasing rcoulomb.
413 set.rcut_coulomb = pme_lb->rcut_coulomb_start;
416 if (pme_lb->cutoff_scheme == CutoffScheme::Verlet)
418 /* Never decrease the Coulomb and VdW list buffers */
419 set.rlistOuter = std::max(set.rcut_coulomb + pme_lb->rbufOuter_coulomb,
420 pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw);
421 set.rlistInner = std::max(set.rcut_coulomb + pme_lb->rbufInner_coulomb,
422 pme_lb->rcut_vdw + pme_lb->rbufInner_vdw);
426 /* TODO Remove these lines and pme_lb->cutoff_scheme */
427 tmpr_coulomb = set.rcut_coulomb + pme_lb->rbufOuter_coulomb;
428 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw;
429 /* Two (known) bugs with cutoff-scheme=group here:
430 * - This modification of rlist results in incorrect DD comunication.
431 * - We should set fr->bTwinRange = (fr->rlistlong > fr->rlist).
433 set.rlistOuter = std::min(tmpr_coulomb, tmpr_vdw);
434 set.rlistInner = set.rlistOuter;
438 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
439 set.grid_efficiency = 1;
440 for (d = 0; d < DIM; d++)
442 set.grid_efficiency *= (set.grid[d] * sp) / norm(pme_lb->box_start[d]);
444 /* The Ewald coefficient is inversly proportional to the cut-off */
445 set.ewaldcoeff_q = pme_lb->setup[0].ewaldcoeff_q * pme_lb->setup[0].rcut_coulomb / set.rcut_coulomb;
446 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
447 set.ewaldcoeff_lj = pme_lb->setup[0].ewaldcoeff_lj * pme_lb->setup[0].rcut_coulomb / set.rcut_coulomb;
455 "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
461 pme_lb->setup.push_back(set);
465 /*! \brief Print the PME grid */
466 static void print_grid(FILE* fp_err, FILE* fp_log, const char* pre, const char* desc, const pme_setup_t* set, double cycles)
468 auto buf = gmx::formatString("%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f",
477 buf += gmx::formatString(": %.1f M-cycles", cycles * 1e-6);
479 if (fp_err != nullptr)
481 fprintf(fp_err, "\r%s\n", buf.c_str());
484 if (fp_log != nullptr)
486 fprintf(fp_log, "%s\n", buf.c_str());
490 /*! \brief Return the index of the last setup used in PME load balancing */
491 static int pme_loadbal_end(pme_load_balancing_t* pme_lb)
493 /* In the initial stage only n is set; end is not set yet */
500 return pme_lb->setup.size();
504 /*! \brief Print descriptive string about what limits PME load balancing */
505 static void print_loadbal_limited(FILE* fp_err, FILE* fp_log, int64_t step, pme_load_balancing_t* pme_lb)
507 auto buf = gmx::formatString(
508 "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
509 gmx::int64ToString(step).c_str(),
510 enumValueToString(pme_lb->elimited),
511 pme_lb->setup[pme_loadbal_end(pme_lb) - 1].rcut_coulomb);
512 if (fp_err != nullptr)
514 fprintf(fp_err, "\r%s\n", buf.c_str());
517 if (fp_log != nullptr)
519 fprintf(fp_log, "%s\n", buf.c_str());
523 /*! \brief Switch load balancing to stage 1
525 * In this stage, only reasonably fast setups are run again. */
526 static void switch_to_stage1(pme_load_balancing_t* pme_lb)
528 /* Increase start until we find a setup that is not slower than
529 * maxRelativeSlowdownAccepted times the fastest setup.
531 pme_lb->start = pme_lb->lower_limit;
532 while (pme_lb->start + 1 < gmx::ssize(pme_lb->setup)
533 && (pme_lb->setup[pme_lb->start].count == 0
534 || pme_lb->setup[pme_lb->start].cycles
535 > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted))
539 /* While increasing start, we might have skipped setups that we did not
540 * time during stage 0. We want to extend the range for stage 1 to include
541 * any skipped setups that lie between setups that were measured to be
542 * acceptably fast and too slow.
544 while (pme_lb->start > pme_lb->lower_limit && pme_lb->setup[pme_lb->start - 1].count == 0)
549 /* Decrease end only with setups that we timed and that are slow. */
550 pme_lb->end = pme_lb->setup.size();
551 if (pme_lb->setup[pme_lb->end - 1].count > 0
552 && pme_lb->setup[pme_lb->end - 1].cycles
553 > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted)
560 /* Next we want to choose setup pme_lb->end-1, but as we will decrease
561 * pme_lb->cur by one right after returning, we set cur to end.
563 pme_lb->cur = pme_lb->end;
566 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
568 * The adjustment is done to generate a different non-bonded PP and PME load.
569 * With separate PME ranks (PP and PME on different processes) or with
570 * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
571 * and changing the load will affect the load balance and performance.
572 * The total time for a set of integration steps is monitored and a range
573 * of grid/cut-off setups is scanned. After calling pme_load_balance many
574 * times and acquiring enough statistics, the best performing setup is chosen.
575 * Here we try to take into account fluctuations and changes due to external
576 * factors as well as DD load balancing.
578 static void pme_load_balance(pme_load_balancing_t* pme_lb,
582 const gmx::MDLogger& mdlog,
583 const t_inputrec& ir,
585 gmx::ArrayRef<const gmx::RVec> x,
587 interaction_const_t* ic,
588 struct nonbonded_verlet_t* nbv,
589 struct gmx_pme_t** pmedata,
595 char buf[STRLEN], sbuf[22];
599 gmx_sumd(1, &cycles, cr);
600 cycles /= cr->nnodes;
603 set = &pme_lb->setup[pme_lb->cur];
606 /* Skip the first c_numPostSwitchTuningIntervalSkip cycles because the first step
607 * after a switch is much slower due to allocation and/or caching effects.
609 if (set->count % (c_numPostSwitchTuningIntervalSkip + 1) != 0)
614 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
615 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
617 GMX_RELEASE_ASSERT(set->count > c_numPostSwitchTuningIntervalSkip, "We should skip cycles");
618 if (set->count == (c_numPostSwitchTuningIntervalSkip + 1))
620 set->cycles = cycles;
624 if (cycles * maxFluctuationAccepted < set->cycles && pme_lb->stage == pme_lb->nstage - 1)
626 /* The performance went up a lot (due to e.g. DD load balancing).
627 * Add a stage, keep the minima, but rescan all setups.
634 "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this "
636 "Increased the number stages to %d"
637 " and ignoring the previous performance\n",
643 maxFluctuationAccepted,
647 set->cycles = std::min(set->cycles, cycles);
650 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
652 pme_lb->fastest = pme_lb->cur;
654 if (haveDDAtomOrdering(*cr))
656 /* We found a new fastest setting, ensure that with subsequent
657 * shorter cut-off's the dynamic load balancing does not make
658 * the use of the current cut-off impossible. This solution is
659 * a trade-off, as the PME load balancing and DD domain size
660 * load balancing can interact in complex ways.
661 * With the Verlet kernels, DD load imbalance will usually be
662 * mainly due to bonded interaction imbalance, which will often
663 * quickly push the domain boundaries beyond the limit for the
664 * optimal, PME load balanced, cut-off. But it could be that
665 * better overal performance can be obtained with a slightly
666 * shorter cut-off and better DD load balancing.
668 set_dd_dlb_max_cutoff(cr, pme_lb->setup[pme_lb->fastest].rlistOuter);
671 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
673 /* Check in stage 0 if we should stop scanning grids.
674 * Stop when the time is more than maxRelativeSlowDownAccepted longer than the fastest.
676 if (pme_lb->stage == 0 && pme_lb->cur > 0
677 && cycles > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted)
679 pme_lb->setup.resize(pme_lb->cur + 1);
680 /* Done with scanning, go to stage 1 */
681 switch_to_stage1(pme_lb);
684 if (pme_lb->stage == 0)
688 gridsize_start = set->grid[XX] * set->grid[YY] * set->grid[ZZ];
692 if (pme_lb->cur + 1 < gmx::ssize(pme_lb->setup))
694 /* We had already generated the next setup */
699 /* Find the next setup */
700 OK = pme_loadbal_increase_cutoff(pme_lb, ir.pme_order, cr->dd);
704 pme_lb->elimited = PmeLoadBalancingLimit::PmeGrid;
709 && pme_lb->setup[pme_lb->cur + 1].spacing > c_maxSpacingScaling * pme_lb->setup[0].spacing)
712 pme_lb->elimited = PmeLoadBalancingLimit::MaxScaling;
715 if (OK && ir.pbcType != PbcType::No)
717 OK = (gmx::square(pme_lb->setup[pme_lb->cur + 1].rlistOuter)
718 <= max_cutoff2(ir.pbcType, box));
721 pme_lb->elimited = PmeLoadBalancingLimit::Box;
729 if (haveDDAtomOrdering(*cr))
731 OK = change_dd_cutoff(cr, box, x, pme_lb->setup[pme_lb->cur].rlistOuter);
734 /* Failed: do not use this setup */
736 pme_lb->elimited = PmeLoadBalancingLimit::DD;
742 /* We hit the upper limit for the cut-off,
743 * the setup should not go further than cur.
745 pme_lb->setup.resize(pme_lb->cur + 1);
746 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
747 /* Switch to the next stage */
748 switch_to_stage1(pme_lb);
751 && !(pme_lb->setup[pme_lb->cur].grid[XX] * pme_lb->setup[pme_lb->cur].grid[YY]
752 * pme_lb->setup[pme_lb->cur].grid[ZZ]
753 < gridsize_start * gridpointsScaleFactor
754 && pme_lb->setup[pme_lb->cur].grid_efficiency
755 < pme_lb->setup[pme_lb->cur - 1].grid_efficiency * relativeEfficiencyFactor));
758 if (pme_lb->stage > 0 && pme_lb->end == 1)
760 pme_lb->cur = pme_lb->lower_limit;
761 pme_lb->stage = pme_lb->nstage;
763 else if (pme_lb->stage > 0 && pme_lb->end > 1)
765 /* If stage = nstage-1:
766 * scan over all setups, rerunning only those setups
767 * which are not much slower than the fastest
770 * Note that we loop backward to minimize the risk of the cut-off
771 * getting limited by DD DLB, since the DLB cut-off limit is set
772 * to the fastest PME setup.
776 if (pme_lb->cur > pme_lb->start)
784 pme_lb->cur = pme_lb->end - 1;
786 } while (pme_lb->stage == pme_lb->nstage - 1 && pme_lb->setup[pme_lb->cur].count > 0
787 && pme_lb->setup[pme_lb->cur].cycles > cycles_fast * maxRelativeSlowdownAccepted);
789 if (pme_lb->stage == pme_lb->nstage)
791 /* We are done optimizing, use the fastest setup we found */
792 pme_lb->cur = pme_lb->fastest;
796 if (haveDDAtomOrdering(*cr) && pme_lb->stage > 0)
798 OK = change_dd_cutoff(cr, box, x, pme_lb->setup[pme_lb->cur].rlistOuter);
801 /* For some reason the chosen cut-off is incompatible with DD.
802 * We should continue scanning a more limited range of cut-off's.
804 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
806 /* stage=nstage says we're finished, but we should continue
807 * balancing, so we set back stage which was just incremented.
811 if (pme_lb->cur <= pme_lb->fastest)
813 /* This should not happen, as we set limits on the DLB bounds.
814 * But we implement a complete failsafe solution anyhow.
816 GMX_LOG(mdlog.warning)
818 .appendTextFormatted(
819 "The fastest PP/PME load balancing setting (cutoff %.3d nm) is no "
820 "longer available due to DD DLB or box size limitations",
822 pme_lb->fastest = pme_lb->lower_limit;
823 pme_lb->start = pme_lb->lower_limit;
825 /* Limit the range to below the current cut-off, scan from start */
826 pme_lb->end = pme_lb->cur;
827 pme_lb->cur = pme_lb->start;
828 pme_lb->elimited = PmeLoadBalancingLimit::DD;
829 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
833 /* Change the Coulomb cut-off and the PME grid */
835 set = &pme_lb->setup[pme_lb->cur];
837 ic->rcoulomb = set->rcut_coulomb;
838 nbv->changePairlistRadii(set->rlistOuter, set->rlistInner);
839 ic->ewaldcoeff_q = set->ewaldcoeff_q;
840 /* TODO: centralize the code that sets the potentials shifts */
841 if (ic->coulomb_modifier == InteractionModifiers::PotShift)
843 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
844 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
846 if (EVDW_PME(ic->vdwtype))
848 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
849 ic->rvdw = set->rcut_coulomb;
850 ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
851 if (ic->vdw_modifier == InteractionModifiers::PotShift)
855 ic->dispersion_shift.cpot = -1.0 / gmx::power6(static_cast<double>(ic->rvdw));
856 ic->repulsion_shift.cpot = -1.0 / gmx::power12(static_cast<double>(ic->rvdw));
857 crc2 = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
859 (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
863 /* We always re-initialize the tables whether they are used or not */
864 init_interaction_const_tables(nullptr, ic, set->rlistOuter, ir.tabext);
866 Nbnxm::gpu_pme_loadbal_update_param(nbv, *ic);
868 if (!pme_lb->bSepPMERanks)
871 * CPU PME keeps a list of allocated pmedata's, that's why pme_lb->setup[pme_lb->cur].pmedata is not always nullptr.
872 * GPU PME, however, currently needs the gmx_pme_reinit always called on load balancing
873 * (pme_gpu_reinit might be not sufficiently decoupled from gmx_pme_init).
874 * This can lead to a lot of reallocations for PME GPU.
875 * Would be nicer if the allocated grid list was hidden within a single pmedata structure.
877 if ((pme_lb->setup[pme_lb->cur].pmedata == nullptr)
878 || pme_gpu_task_enabled(pme_lb->setup[pme_lb->cur].pmedata))
880 /* Generate a new PME data structure,
881 * copying part of the old pointers.
884 &set->pmedata, cr, pme_lb->setup[0].pmedata, &ir, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
886 *pmedata = set->pmedata;
890 /* Tell our PME-only rank to switch grid */
891 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
896 print_grid(nullptr, debug, "", "switched to", set, -1);
899 if (pme_lb->stage == pme_lb->nstage)
901 print_grid(fp_err, fp_log, "", "optimal", set, -1);
905 /*! \brief Prepare for another round of PME load balancing
907 * \param[in,out] pme_lb Pointer to PME load balancing struct
908 * \param[in] bDlbUnlocked TRUE is DLB was locked and is now unlocked
910 * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
911 * the PP/PME balance might change and re-balancing can improve performance.
912 * This function adds 2 stages and adjusts the considered setup range.
914 static void continue_pme_loadbal(pme_load_balancing_t* pme_lb, gmx_bool bDlbUnlocked)
916 /* Add 2 tuning stages, keep the detected end of the setup range */
918 if (bDlbUnlocked && pme_lb->bSepPMERanks)
920 /* With separate PME ranks, DLB should always lower the PP load and
921 * can only increase the PME load (more communication and imbalance),
922 * so we only need to scan longer cut-off's.
924 pme_lb->lower_limit = pme_lb->cur;
926 pme_lb->start = pme_lb->lower_limit;
929 void pme_loadbal_do(pme_load_balancing_t* pme_lb,
933 const gmx::MDLogger& mdlog,
934 const t_inputrec& ir,
937 gmx::ArrayRef<const gmx::RVec> x,
938 gmx_wallcycle* wcycle,
942 bool useGpuPmePpCommunication)
947 assert(pme_lb != nullptr);
949 if (!pme_lb->bActive)
954 n_prev = pme_lb->cycles_n;
955 cycles_prev = pme_lb->cycles_c;
956 wallcycle_get(wcycle, WallCycleCounter::Step, &pme_lb->cycles_n, &pme_lb->cycles_c);
958 /* Before the first step we haven't done any steps yet.
959 * Also handle cases where ir.init_step % ir.nstlist != 0.
960 * We also want to skip a number of steps and seconds while
961 * the CPU and GPU, when used, performance stabilizes.
963 if (!PAR(cr) || (haveDDAtomOrdering(*cr) && DDMASTER(cr->dd)))
965 pme_lb->startupTimeDelayElapsed = (gmx_gettime() - pme_lb->startTime < c_startupTimeDelay);
967 if (haveDDAtomOrdering(*cr))
969 dd_bcast(cr->dd, sizeof(bool), &pme_lb->startupTimeDelayElapsed);
972 if (pme_lb->cycles_n == 0 || step_rel < c_numFirstTuningIntervalSkip * ir.nstlist
973 || pme_lb->startupTimeDelayElapsed)
978 /* Sanity check, we expect nstlist cycle counts */
979 if (pme_lb->cycles_n - n_prev != ir.nstlist)
981 /* We could return here, but it's safer to issue an error and quit */
982 gmx_incons("pme_loadbal_do called at an interval != nstlist");
985 /* PME grid + cut-off optimization with GPUs or PME ranks */
986 if (!pme_lb->bBalance && pme_lb->bSepPMERanks)
988 if (pme_lb->bTriggerOnDLB)
990 pme_lb->bBalance = dd_dlb_is_on(cr->dd);
992 /* We should ignore the first timing to avoid timing allocation
993 * overhead. And since the PME load balancing is called just
994 * before DD repartitioning, the ratio returned by dd_pme_f_ratio
995 * is not over the last nstlist steps, but the nstlist steps before
996 * that. So the first useful ratio is available at step_rel=3*nstlist.
998 else if (step_rel >= c_numFirstTuningIntervalSkipWithSepPme * ir.nstlist)
1000 GMX_ASSERT(haveDDAtomOrdering(*cr), "Domain decomposition should be active here");
1001 if (DDMASTER(cr->dd))
1003 /* If PME rank load is too high, start tuning. If
1004 PME-PP direct GPU communication is active,
1005 unconditionally start tuning since ratio will be
1006 unreliable due to CPU-GPU asynchronicity in codepath */
1007 pme_lb->bBalance = useGpuPmePpCommunication
1009 : (dd_pme_f_ratio(cr->dd) >= loadBalanceTriggerFactor);
1011 dd_bcast(cr->dd, sizeof(gmx_bool), &pme_lb->bBalance);
1014 pme_lb->bActive = (pme_lb->bBalance || step_rel <= pme_lb->step_rel_stop);
1017 /* The location in the code of this balancing termination is strange.
1018 * You would expect to have it after the call to pme_load_balance()
1019 * below, since there pme_lb->stage is updated.
1020 * But when terminating directly after deciding on and selecting the
1021 * optimal setup, DLB will turn on right away if it was locked before.
1022 * This might be due to PME reinitialization. So we check stage here
1023 * to allow for another nstlist steps with DLB locked to stabilize
1026 if (pme_lb->bBalance && pme_lb->stage == pme_lb->nstage)
1028 pme_lb->bBalance = FALSE;
1030 if (haveDDAtomOrdering(*cr) && dd_dlb_is_locked(cr->dd))
1032 /* Unlock the DLB=auto, DLB is allowed to activate */
1033 dd_dlb_unlock(cr->dd);
1034 GMX_LOG(mdlog.warning)
1036 .appendText("NOTE: DLB can now turn on, when beneficial");
1038 /* We don't deactivate the tuning yet, since we will balance again
1039 * after DLB gets turned on, if it does within PMETune_period.
1041 continue_pme_loadbal(pme_lb, TRUE);
1042 pme_lb->bTriggerOnDLB = TRUE;
1043 pme_lb->step_rel_stop = step_rel + PMETunePeriod * ir.nstlist;
1047 /* We're completely done with PME tuning */
1048 pme_lb->bActive = FALSE;
1051 if (haveDDAtomOrdering(*cr))
1053 /* Set the cut-off limit to the final selected cut-off,
1054 * so we don't have artificial DLB limits.
1055 * This also ensures that we won't disable the currently
1056 * optimal setting during a second round of PME balancing.
1058 set_dd_dlb_max_cutoff(cr, fr->nbv->pairlistOuterRadius());
1062 if (pme_lb->bBalance)
1064 /* We might not have collected nstlist steps in cycles yet,
1065 * since init_step might not be a multiple of nstlist,
1066 * but the first data collected is skipped anyhow.
1068 pme_load_balance(pme_lb,
1076 pme_lb->cycles_c - cycles_prev,
1082 /* Update deprecated rlist in forcerec to stay in sync with fr->nbv */
1083 fr->rlist = fr->nbv->pairlistOuterRadius();
1085 if (ir.eDispCorr != DispersionCorrectionType::No)
1087 fr->dispersionCorrection->setParameters(*fr->ic);
1091 if (!pme_lb->bBalance && (!pme_lb->bSepPMERanks || step_rel > pme_lb->step_rel_stop))
1093 /* We have just deactivated the balancing and we're not measuring PP/PME
1094 * imbalance during the first steps of the run: deactivate the tuning.
1096 pme_lb->bActive = FALSE;
1099 if (!(pme_lb->bActive) && haveDDAtomOrdering(*cr) && dd_dlb_is_locked(cr->dd))
1101 /* Make sure DLB is allowed when we deactivate PME tuning */
1102 dd_dlb_unlock(cr->dd);
1103 GMX_LOG(mdlog.warning)
1105 .appendText("NOTE: DLB can now turn on, when beneficial");
1108 *bPrinting = pme_lb->bBalance;
1111 /*! \brief Return product of the number of PME grid points in each dimension */
1112 static int pme_grid_points(const pme_setup_t* setup)
1114 return setup->grid[XX] * setup->grid[YY] * setup->grid[ZZ];
1117 /*! \brief Print one load-balancing setting */
1118 static void print_pme_loadbal_setting(FILE* fplog, const char* name, const pme_setup_t* setup)
1121 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
1123 setup->rcut_coulomb,
1129 1 / setup->ewaldcoeff_q);
1132 /*! \brief Print all load-balancing settings */
1133 static void print_pme_loadbal_settings(pme_load_balancing_t* pme_lb,
1135 const gmx::MDLogger& mdlog,
1136 gmx_bool bNonBondedOnGPU)
1138 double pp_ratio, grid_ratio;
1139 real pp_ratio_temporary;
1141 pp_ratio_temporary = pme_lb->setup[pme_lb->cur].rlistInner / pme_lb->setup[0].rlistInner;
1142 pp_ratio = gmx::power3(pp_ratio_temporary);
1143 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])
1144 / static_cast<double>(pme_grid_points(&pme_lb->setup[0]));
1146 fprintf(fplog, "\n");
1147 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
1148 fprintf(fplog, "\n");
1149 /* Here we only warn when the optimal setting is the last one */
1150 if (pme_lb->elimited != PmeLoadBalancingLimit::No && pme_lb->cur == pme_loadbal_end(pme_lb) - 1)
1153 " NOTE: The PP/PME load balancing was limited by the %s,\n",
1154 enumValueToString(pme_lb->elimited));
1155 fprintf(fplog, " you might not have reached a good load balance.\n");
1156 if (pme_lb->elimited == PmeLoadBalancingLimit::DD)
1158 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
1160 fprintf(fplog, "\n");
1162 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
1163 fprintf(fplog, " particle-particle PME\n");
1164 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
1165 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
1166 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
1167 fprintf(fplog, " cost-ratio %4.2f %4.2f\n", pp_ratio, grid_ratio);
1168 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
1170 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
1172 GMX_LOG(mdlog.warning)
1175 "NOTE: PME load balancing increased the non-bonded workload by more than "
1177 " For better performance, use (more) PME ranks (mdrun -npme),\n"
1178 " or if you are beyond the scaling limit, use fewer total ranks (or "
1183 fprintf(fplog, "\n");
1187 void pme_loadbal_done(pme_load_balancing_t* pme_lb, FILE* fplog, const gmx::MDLogger& mdlog, gmx_bool bNonBondedOnGPU)
1189 if (fplog != nullptr && (pme_lb->cur > 0 || pme_lb->elimited != PmeLoadBalancingLimit::No))
1191 print_pme_loadbal_settings(pme_lb, fplog, mdlog, bNonBondedOnGPU);