2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5 * Copyright (c) 2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
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 "pme_load_balancing.h"
54 #include "gromacs/domdec/dlb.h"
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_network.h"
57 #include "gromacs/domdec/domdec_struct.h"
58 #include "gromacs/domdec/partition.h"
59 #include "gromacs/ewald/ewald_utils.h"
60 #include "gromacs/ewald/pme.h"
61 #include "gromacs/fft/calcgrid.h"
62 #include "gromacs/gmxlib/network.h"
63 #include "gromacs/math/functions.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/mdlib/dispersioncorrection.h"
66 #include "gromacs/mdlib/forcerec.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/nbnxm/gpu_data_mgmt.h"
72 #include "gromacs/nbnxm/nbnxm.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/logger.h"
79 #include "gromacs/utility/smalloc.h"
80 #include "gromacs/utility/strconvert.h"
82 #include "pme_internal.h"
84 /*! \brief Parameters and settings for one PP-PME setup */
87 real rcut_coulomb; /**< Coulomb cut-off */
88 real rlistOuter; /**< cut-off for the outer pair-list */
89 real rlistInner; /**< cut-off for the inner pair-list */
90 real spacing; /**< (largest) PME grid spacing */
91 ivec grid; /**< the PME grid dimensions */
92 real grid_efficiency; /**< ineffiency factor for non-uniform grids <= 1 */
93 real ewaldcoeff_q; /**< Electrostatic Ewald coefficient */
94 real ewaldcoeff_lj; /**< LJ Ewald coefficient, only for the call to send_switchgrid */
95 struct gmx_pme_t* pmedata; /**< the data structure used in the PME code */
96 int count; /**< number of times this setup has been timed */
97 double cycles; /**< the fastest time for this setup in cycles */
100 /*! \brief After 50 nstlist periods of not observing imbalance: never tune PME */
101 const int PMETunePeriod = 50;
102 /*! \brief Trigger PME load balancing at more than 5% PME overload */
103 const real loadBalanceTriggerFactor = 1.05;
104 /*! \brief Scale the grid by a most at factor 1.7.
106 * This still leaves room for about 4-4.5x decrease in grid spacing while limiting the cases where
107 * large imbalance leads to extreme cutoff scaling for marginal benefits.
109 * This should help to avoid:
110 * - large increase in power consumption for little performance gain
111 * - increasing communication volume
114 const real c_maxSpacingScaling = 1.7;
115 /*! \brief In the initial scan, step by grids that are at least a factor 0.8 coarser */
116 const real gridpointsScaleFactor = 0.8;
117 /*! \brief In the initial scan, try to skip grids with uneven x/y/z spacing,
118 * checking if the "efficiency" is more than 5% worse than the previous grid.
120 const real relativeEfficiencyFactor = 1.05;
121 /*! \brief Rerun until a run is 12% slower setups than the fastest run so far */
122 const real maxRelativeSlowdownAccepted = 1.12;
123 /*! \brief If setups get more than 2% faster, do another round to avoid
124 * choosing a slower setup due to acceleration or fluctuations.
126 const real maxFluctuationAccepted = 1.02;
128 //! \brief Number of nstlist long tuning intervals to skip before starting
129 // load-balancing at the beginning of the run.
130 const int c_numFirstTuningIntervalSkip = 5;
131 //! \brief Number of nstlist long tuning intervals to skip before starting
132 // load-balancing at the beginning of the run with separate PME ranks. */
133 const int c_numFirstTuningIntervalSkipWithSepPme = 3;
134 //! \brief Number of nstlist long tuning intervals to skip after switching to a new setting
136 const int c_numPostSwitchTuningIntervalSkip = 1;
137 //! \brief Number of seconds to delay the tuning at startup to allow processors clocks to ramp up.
138 const double c_startupTimeDelay = 5.0;
140 /*! \brief Enumeration whose values describe the effect limiting the load balancing */
151 /*! \brief Descriptive strings matching ::epmelb */
152 static const char* pmelblim_str[epmelblimNR] = { "no", "box size", "domain decompostion",
153 "PME grid restriction",
154 "maximum allowed grid scaling" };
156 struct pme_load_balancing_t
158 gmx_bool bSepPMERanks; /**< do we have separate PME ranks? */
159 gmx_bool bActive; /**< is PME tuning active? */
160 int64_t step_rel_stop; /**< stop the tuning after this value of step_rel */
161 gmx_bool bTriggerOnDLB; /**< trigger balancing only on DD DLB */
162 gmx_bool bBalance; /**< are we in the balancing phase, i.e. trying different setups? */
163 int nstage; /**< the current maximum number of stages */
164 bool startupTimeDelayElapsed; /**< Has the c_startupTimeDelay elapsed indicating that the balancing can start. */
166 real cut_spacing; /**< the minimum cutoff / PME grid spacing ratio */
167 real rcut_vdw; /**< Vdw cutoff (does not change) */
168 real rcut_coulomb_start; /**< Initial electrostatics cutoff */
169 real rbufOuter_coulomb; /**< the outer pairlist buffer size */
170 real rbufOuter_vdw; /**< the outer pairlist buffer size */
171 real rbufInner_coulomb; /**< the inner pairlist buffer size */
172 real rbufInner_vdw; /**< the inner pairlist buffer size */
173 matrix box_start; /**< the initial simulation box */
174 std::vector<pme_setup_t> setup; /**< the PME+cutoff setups */
175 int cur; /**< the index (in setup) of the current setup */
176 int fastest; /**< index of the fastest setup up till now */
177 int lower_limit; /**< don't go below this setup index */
178 int start; /**< start of setup index range to consider in stage>0 */
179 int end; /**< end of setup index range to consider in stage>0 */
180 int elimited; /**< was the balancing limited, uses enum above */
181 int cutoff_scheme; /**< Verlet or group cut-offs */
183 int stage; /**< the current stage */
185 int cycles_n; /**< step cycle counter cumulative count */
186 double cycles_c; /**< step cycle counter cumulative cycles */
187 double startTime; /**< time stamp when the balancing was started on the master rank (relative to the UNIX epoch start).*/
190 /* TODO The code in this file should call this getter, rather than
191 * read bActive anywhere */
192 bool pme_loadbal_is_active(const pme_load_balancing_t* pme_lb)
194 return pme_lb != nullptr && pme_lb->bActive;
197 // TODO Return a unique_ptr to pme_load_balancing_t
198 void pme_loadbal_init(pme_load_balancing_t** pme_lb_p,
200 const gmx::MDLogger& mdlog,
201 const t_inputrec& ir,
203 const interaction_const_t& ic,
204 const nonbonded_verlet_t& nbv,
209 pme_load_balancing_t* pme_lb;
213 // Note that we don't (yet) support PME load balancing with LJ-PME only.
214 GMX_RELEASE_ASSERT(EEL_PME(ir.coulombtype),
215 "pme_loadbal_init called without PME electrostatics");
216 // To avoid complexity, we require a single cut-off with PME for q+LJ.
217 // This is checked by grompp, but it doesn't hurt to check again.
218 GMX_RELEASE_ASSERT(!(EEL_PME(ir.coulombtype) && EVDW_PME(ir.vdwtype) && ir.rcoulomb != ir.rvdw),
219 "With Coulomb and LJ PME, rcoulomb should be equal to rvdw");
221 pme_lb = new pme_load_balancing_t;
223 pme_lb->bSepPMERanks = !thisRankHasDuty(cr, DUTY_PME);
225 /* Initially we turn on balancing directly on based on PP/PME imbalance */
226 pme_lb->bTriggerOnDLB = FALSE;
228 /* Any number of stages >= 2 is supported */
231 pme_lb->cutoff_scheme = ir.cutoff_scheme;
233 pme_lb->rbufOuter_coulomb = nbv.pairlistOuterRadius() - ic.rcoulomb;
234 pme_lb->rbufOuter_vdw = nbv.pairlistOuterRadius() - ic.rvdw;
235 pme_lb->rbufInner_coulomb = nbv.pairlistInnerRadius() - ic.rcoulomb;
236 pme_lb->rbufInner_vdw = nbv.pairlistInnerRadius() - ic.rvdw;
238 /* Scale box with Ewald wall factor; note that we pmedata->boxScaler
239 * can't always usedd as it's not available with separate PME ranks.
241 EwaldBoxZScaler boxScaler(ir);
242 boxScaler.scaleBox(box, pme_lb->box_start);
244 pme_lb->setup.resize(1);
246 pme_lb->rcut_vdw = ic.rvdw;
247 pme_lb->rcut_coulomb_start = ir.rcoulomb;
250 pme_lb->setup[0].rcut_coulomb = ic.rcoulomb;
251 pme_lb->setup[0].rlistOuter = nbv.pairlistOuterRadius();
252 pme_lb->setup[0].rlistInner = nbv.pairlistInnerRadius();
253 pme_lb->setup[0].grid[XX] = ir.nkx;
254 pme_lb->setup[0].grid[YY] = ir.nky;
255 pme_lb->setup[0].grid[ZZ] = ir.nkz;
256 pme_lb->setup[0].ewaldcoeff_q = ic.ewaldcoeff_q;
257 pme_lb->setup[0].ewaldcoeff_lj = ic.ewaldcoeff_lj;
259 if (!pme_lb->bSepPMERanks)
261 GMX_RELEASE_ASSERT(pmedata,
262 "On ranks doing both PP and PME we need a valid pmedata object");
263 pme_lb->setup[0].pmedata = pmedata;
267 for (d = 0; d < DIM; d++)
269 sp = norm(pme_lb->box_start[d]) / pme_lb->setup[0].grid[d];
275 pme_lb->setup[0].spacing = spm;
277 if (ir.fourier_spacing > 0)
279 pme_lb->cut_spacing = ir.rcoulomb / ir.fourier_spacing;
283 pme_lb->cut_spacing = ir.rcoulomb / pme_lb->setup[0].spacing;
289 pme_lb->lower_limit = 0;
292 pme_lb->elimited = epmelblimNO;
294 pme_lb->cycles_n = 0;
295 pme_lb->cycles_c = 0;
296 // only master ranks do timing
297 if (!PAR(cr) || (DOMAINDECOMP(cr) && DDMASTER(cr->dd)))
299 pme_lb->startTime = gmx_gettime();
302 if (!wallcycle_have_counter())
304 GMX_LOG(mdlog.warning)
307 "NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use "
308 "PME-PP balancing.");
311 /* Tune with GPUs and/or separate PME ranks.
312 * When running only on a CPU without PME ranks, PME tuning will only help
313 * with small numbers of atoms in the cut-off sphere.
315 pme_lb->bActive = (wallcycle_have_counter() && (bUseGPU || pme_lb->bSepPMERanks));
317 /* With GPUs and no separate PME ranks we can't measure the PP/PME
318 * imbalance, so we start balancing right away.
319 * Otherwise we only start balancing after we observe imbalance.
321 pme_lb->bBalance = (pme_lb->bActive && (bUseGPU && !pme_lb->bSepPMERanks));
323 pme_lb->step_rel_stop = PMETunePeriod * ir.nstlist;
325 /* Delay DD load balancing when GPUs are used */
326 if (pme_lb->bActive && DOMAINDECOMP(cr) && cr->dd->nnodes > 1 && bUseGPU)
328 /* Lock DLB=auto to off (does nothing when DLB=yes/no.
329 * With GPUs and separate PME nodes, we want to first
330 * do PME tuning without DLB, since DLB might limit
331 * the cut-off, which never improves performance.
332 * We allow for DLB + PME tuning after a first round of tuning.
335 if (dd_dlb_is_locked(cr->dd))
337 GMX_LOG(mdlog.warning)
339 .appendText("NOTE: DLB will not turn on during the first phase of PME tuning");
346 /*! \brief Try to increase the cutoff during load balancing */
347 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t* pme_lb, int pme_order, const gmx_domdec_t* dd)
350 real tmpr_coulomb, tmpr_vdw;
354 /* Try to add a new setup with next larger cut-off to the list */
357 set.pmedata = nullptr;
359 NumPmeDomains numPmeDomains = getNumPmeDomains(dd);
364 /* Avoid infinite while loop, which can occur at the minimum grid size.
365 * Note that in practice load balancing will stop before this point.
366 * The factor 2.1 allows for the extreme case in which only grids
367 * of powers of 2 are allowed (the current code supports more grids).
375 clear_ivec(set.grid);
376 sp = calcFftGrid(nullptr, pme_lb->box_start, fac * pme_lb->setup[pme_lb->cur].spacing,
377 minimalPmeGridSize(pme_order), &set.grid[XX], &set.grid[YY], &set.grid[ZZ]);
379 /* As here we can't easily check if one of the PME ranks
380 * uses threading, we do a conservative grid check.
381 * This means we can't use pme_order or less grid lines
382 * per PME rank along x, which is not a strong restriction.
384 grid_ok = gmx_pme_check_restrictions(pme_order, set.grid[XX], set.grid[YY], set.grid[ZZ],
385 numPmeDomains.x, true, false);
386 } while (sp <= 1.001 * pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
388 set.rcut_coulomb = pme_lb->cut_spacing * sp;
389 if (set.rcut_coulomb < pme_lb->rcut_coulomb_start)
391 /* This is unlikely, but can happen when e.g. continuing from
392 * a checkpoint after equilibration where the box shrank a lot.
393 * We want to avoid rcoulomb getting smaller than rvdw
394 * and there might be more issues with decreasing rcoulomb.
396 set.rcut_coulomb = pme_lb->rcut_coulomb_start;
399 if (pme_lb->cutoff_scheme == ecutsVERLET)
401 /* Never decrease the Coulomb and VdW list buffers */
402 set.rlistOuter = std::max(set.rcut_coulomb + pme_lb->rbufOuter_coulomb,
403 pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw);
404 set.rlistInner = std::max(set.rcut_coulomb + pme_lb->rbufInner_coulomb,
405 pme_lb->rcut_vdw + pme_lb->rbufInner_vdw);
409 /* TODO Remove these lines and pme_lb->cutoff_scheme */
410 tmpr_coulomb = set.rcut_coulomb + pme_lb->rbufOuter_coulomb;
411 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw;
412 /* Two (known) bugs with cutoff-scheme=group here:
413 * - This modification of rlist results in incorrect DD comunication.
414 * - We should set fr->bTwinRange = (fr->rlistlong > fr->rlist).
416 set.rlistOuter = std::min(tmpr_coulomb, tmpr_vdw);
417 set.rlistInner = set.rlistOuter;
421 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
422 set.grid_efficiency = 1;
423 for (d = 0; d < DIM; d++)
425 set.grid_efficiency *= (set.grid[d] * sp) / norm(pme_lb->box_start[d]);
427 /* The Ewald coefficient is inversly proportional to the cut-off */
428 set.ewaldcoeff_q = pme_lb->setup[0].ewaldcoeff_q * pme_lb->setup[0].rcut_coulomb / set.rcut_coulomb;
429 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
430 set.ewaldcoeff_lj = pme_lb->setup[0].ewaldcoeff_lj * pme_lb->setup[0].rcut_coulomb / set.rcut_coulomb;
437 fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n", set.grid[XX],
438 set.grid[YY], set.grid[ZZ], set.rcut_coulomb);
440 pme_lb->setup.push_back(set);
444 /*! \brief Print the PME grid */
445 static void print_grid(FILE* fp_err, FILE* fp_log, const char* pre, const char* desc, const pme_setup_t* set, double cycles)
447 auto buf = gmx::formatString("%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f", pre, desc,
448 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
451 buf += gmx::formatString(": %.1f M-cycles", cycles * 1e-6);
453 if (fp_err != nullptr)
455 fprintf(fp_err, "\r%s\n", buf.c_str());
458 if (fp_log != nullptr)
460 fprintf(fp_log, "%s\n", buf.c_str());
464 /*! \brief Return the index of the last setup used in PME load balancing */
465 static int pme_loadbal_end(pme_load_balancing_t* pme_lb)
467 /* In the initial stage only n is set; end is not set yet */
474 return pme_lb->setup.size();
478 /*! \brief Print descriptive string about what limits PME load balancing */
479 static void print_loadbal_limited(FILE* fp_err, FILE* fp_log, int64_t step, pme_load_balancing_t* pme_lb)
481 auto buf = gmx::formatString(
482 "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
483 gmx::int64ToString(step).c_str(), pmelblim_str[pme_lb->elimited],
484 pme_lb->setup[pme_loadbal_end(pme_lb) - 1].rcut_coulomb);
485 if (fp_err != nullptr)
487 fprintf(fp_err, "\r%s\n", buf.c_str());
490 if (fp_log != nullptr)
492 fprintf(fp_log, "%s\n", buf.c_str());
496 /*! \brief Switch load balancing to stage 1
498 * In this stage, only reasonably fast setups are run again. */
499 static void switch_to_stage1(pme_load_balancing_t* pme_lb)
501 /* Increase start until we find a setup that is not slower than
502 * maxRelativeSlowdownAccepted times the fastest setup.
504 pme_lb->start = pme_lb->lower_limit;
505 while (pme_lb->start + 1 < gmx::ssize(pme_lb->setup)
506 && (pme_lb->setup[pme_lb->start].count == 0
507 || pme_lb->setup[pme_lb->start].cycles
508 > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted))
512 /* While increasing start, we might have skipped setups that we did not
513 * time during stage 0. We want to extend the range for stage 1 to include
514 * any skipped setups that lie between setups that were measured to be
515 * acceptably fast and too slow.
517 while (pme_lb->start > pme_lb->lower_limit && pme_lb->setup[pme_lb->start - 1].count == 0)
522 /* Decrease end only with setups that we timed and that are slow. */
523 pme_lb->end = pme_lb->setup.size();
524 if (pme_lb->setup[pme_lb->end - 1].count > 0
525 && pme_lb->setup[pme_lb->end - 1].cycles
526 > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted)
533 /* Next we want to choose setup pme_lb->end-1, but as we will decrease
534 * pme_lb->cur by one right after returning, we set cur to end.
536 pme_lb->cur = pme_lb->end;
539 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
541 * The adjustment is done to generate a different non-bonded PP and PME load.
542 * With separate PME ranks (PP and PME on different processes) or with
543 * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
544 * and changing the load will affect the load balance and performance.
545 * The total time for a set of integration steps is monitored and a range
546 * of grid/cut-off setups is scanned. After calling pme_load_balance many
547 * times and acquiring enough statistics, the best performing setup is chosen.
548 * Here we try to take into account fluctuations and changes due to external
549 * factors as well as DD load balancing.
551 static void pme_load_balance(pme_load_balancing_t* pme_lb,
555 const gmx::MDLogger& mdlog,
556 const t_inputrec& ir,
558 gmx::ArrayRef<const gmx::RVec> x,
560 interaction_const_t* ic,
561 struct nonbonded_verlet_t* nbv,
562 struct gmx_pme_t** pmedata,
568 char buf[STRLEN], sbuf[22];
572 gmx_sumd(1, &cycles, cr);
573 cycles /= cr->nnodes;
576 set = &pme_lb->setup[pme_lb->cur];
579 /* Skip the first c_numPostSwitchTuningIntervalSkip cycles because the first step
580 * after a switch is much slower due to allocation and/or caching effects.
582 if (set->count % (c_numPostSwitchTuningIntervalSkip + 1) != 0)
587 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
588 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
590 GMX_RELEASE_ASSERT(set->count > c_numPostSwitchTuningIntervalSkip, "We should skip cycles");
591 if (set->count == (c_numPostSwitchTuningIntervalSkip + 1))
593 set->cycles = cycles;
597 if (cycles * maxFluctuationAccepted < set->cycles && pme_lb->stage == pme_lb->nstage - 1)
599 /* The performance went up a lot (due to e.g. DD load balancing).
600 * Add a stage, keep the minima, but rescan all setups.
607 "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this "
609 "Increased the number stages to %d"
610 " and ignoring the previous performance\n",
611 set->grid[XX], set->grid[YY], set->grid[ZZ], set->cycles * 1e-6,
612 cycles * 1e-6, maxFluctuationAccepted, pme_lb->nstage);
615 set->cycles = std::min(set->cycles, cycles);
618 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
620 pme_lb->fastest = pme_lb->cur;
622 if (DOMAINDECOMP(cr))
624 /* We found a new fastest setting, ensure that with subsequent
625 * shorter cut-off's the dynamic load balancing does not make
626 * the use of the current cut-off impossible. This solution is
627 * a trade-off, as the PME load balancing and DD domain size
628 * load balancing can interact in complex ways.
629 * With the Verlet kernels, DD load imbalance will usually be
630 * mainly due to bonded interaction imbalance, which will often
631 * quickly push the domain boundaries beyond the limit for the
632 * optimal, PME load balanced, cut-off. But it could be that
633 * better overal performance can be obtained with a slightly
634 * shorter cut-off and better DD load balancing.
636 set_dd_dlb_max_cutoff(cr, pme_lb->setup[pme_lb->fastest].rlistOuter);
639 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
641 /* Check in stage 0 if we should stop scanning grids.
642 * Stop when the time is more than maxRelativeSlowDownAccepted longer than the fastest.
644 if (pme_lb->stage == 0 && pme_lb->cur > 0
645 && cycles > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted)
647 pme_lb->setup.resize(pme_lb->cur + 1);
648 /* Done with scanning, go to stage 1 */
649 switch_to_stage1(pme_lb);
652 if (pme_lb->stage == 0)
656 gridsize_start = set->grid[XX] * set->grid[YY] * set->grid[ZZ];
660 if (pme_lb->cur + 1 < gmx::ssize(pme_lb->setup))
662 /* We had already generated the next setup */
667 /* Find the next setup */
668 OK = pme_loadbal_increase_cutoff(pme_lb, ir.pme_order, cr->dd);
672 pme_lb->elimited = epmelblimPMEGRID;
677 && pme_lb->setup[pme_lb->cur + 1].spacing > c_maxSpacingScaling * pme_lb->setup[0].spacing)
680 pme_lb->elimited = epmelblimMAXSCALING;
683 if (OK && ir.ePBC != epbcNONE)
685 OK = (gmx::square(pme_lb->setup[pme_lb->cur + 1].rlistOuter) <= max_cutoff2(ir.ePBC, box));
688 pme_lb->elimited = epmelblimBOX;
696 if (DOMAINDECOMP(cr))
698 OK = change_dd_cutoff(cr, box, x, pme_lb->setup[pme_lb->cur].rlistOuter);
701 /* Failed: do not use this setup */
703 pme_lb->elimited = epmelblimDD;
709 /* We hit the upper limit for the cut-off,
710 * the setup should not go further than cur.
712 pme_lb->setup.resize(pme_lb->cur + 1);
713 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
714 /* Switch to the next stage */
715 switch_to_stage1(pme_lb);
718 && !(pme_lb->setup[pme_lb->cur].grid[XX] * pme_lb->setup[pme_lb->cur].grid[YY]
719 * pme_lb->setup[pme_lb->cur].grid[ZZ]
720 < gridsize_start * gridpointsScaleFactor
721 && pme_lb->setup[pme_lb->cur].grid_efficiency
722 < pme_lb->setup[pme_lb->cur - 1].grid_efficiency * relativeEfficiencyFactor));
725 if (pme_lb->stage > 0 && pme_lb->end == 1)
727 pme_lb->cur = pme_lb->lower_limit;
728 pme_lb->stage = pme_lb->nstage;
730 else if (pme_lb->stage > 0 && pme_lb->end > 1)
732 /* If stage = nstage-1:
733 * scan over all setups, rerunning only those setups
734 * which are not much slower than the fastest
737 * Note that we loop backward to minimize the risk of the cut-off
738 * getting limited by DD DLB, since the DLB cut-off limit is set
739 * to the fastest PME setup.
743 if (pme_lb->cur > pme_lb->start)
751 pme_lb->cur = pme_lb->end - 1;
753 } while (pme_lb->stage == pme_lb->nstage - 1 && pme_lb->setup[pme_lb->cur].count > 0
754 && pme_lb->setup[pme_lb->cur].cycles > cycles_fast * maxRelativeSlowdownAccepted);
756 if (pme_lb->stage == pme_lb->nstage)
758 /* We are done optimizing, use the fastest setup we found */
759 pme_lb->cur = pme_lb->fastest;
763 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
765 OK = change_dd_cutoff(cr, box, x, pme_lb->setup[pme_lb->cur].rlistOuter);
768 /* For some reason the chosen cut-off is incompatible with DD.
769 * We should continue scanning a more limited range of cut-off's.
771 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
773 /* stage=nstage says we're finished, but we should continue
774 * balancing, so we set back stage which was just incremented.
778 if (pme_lb->cur <= pme_lb->fastest)
780 /* This should not happen, as we set limits on the DLB bounds.
781 * But we implement a complete failsafe solution anyhow.
783 GMX_LOG(mdlog.warning)
785 .appendTextFormatted(
786 "The fastest PP/PME load balancing setting (cutoff %.3d nm) is no "
787 "longer available due to DD DLB or box size limitations",
789 pme_lb->fastest = pme_lb->lower_limit;
790 pme_lb->start = pme_lb->lower_limit;
792 /* Limit the range to below the current cut-off, scan from start */
793 pme_lb->end = pme_lb->cur;
794 pme_lb->cur = pme_lb->start;
795 pme_lb->elimited = epmelblimDD;
796 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
800 /* Change the Coulomb cut-off and the PME grid */
802 set = &pme_lb->setup[pme_lb->cur];
804 ic->rcoulomb = set->rcut_coulomb;
805 nbv->changePairlistRadii(set->rlistOuter, set->rlistInner);
806 ic->ewaldcoeff_q = set->ewaldcoeff_q;
807 /* TODO: centralize the code that sets the potentials shifts */
808 if (ic->coulomb_modifier == eintmodPOTSHIFT)
810 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
811 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
813 if (EVDW_PME(ic->vdwtype))
815 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
816 ic->rvdw = set->rcut_coulomb;
817 ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
818 if (ic->vdw_modifier == eintmodPOTSHIFT)
822 ic->dispersion_shift.cpot = -1.0 / gmx::power6(static_cast<double>(ic->rvdw));
823 ic->repulsion_shift.cpot = -1.0 / gmx::power12(static_cast<double>(ic->rvdw));
824 crc2 = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
826 (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
830 /* We always re-initialize the tables whether they are used or not */
831 init_interaction_const_tables(nullptr, ic);
833 Nbnxm::gpu_pme_loadbal_update_param(nbv, ic);
835 if (!pme_lb->bSepPMERanks)
838 * CPU PME keeps a list of allocated pmedata's, that's why pme_lb->setup[pme_lb->cur].pmedata is not always nullptr.
839 * GPU PME, however, currently needs the gmx_pme_reinit always called on load balancing
840 * (pme_gpu_reinit might be not sufficiently decoupled from gmx_pme_init).
841 * This can lead to a lot of reallocations for PME GPU.
842 * Would be nicer if the allocated grid list was hidden within a single pmedata structure.
844 if ((pme_lb->setup[pme_lb->cur].pmedata == nullptr)
845 || pme_gpu_task_enabled(pme_lb->setup[pme_lb->cur].pmedata))
847 /* Generate a new PME data structure,
848 * copying part of the old pointers.
850 gmx_pme_reinit(&set->pmedata, cr, pme_lb->setup[0].pmedata, &ir, set->grid,
851 set->ewaldcoeff_q, set->ewaldcoeff_lj);
853 *pmedata = set->pmedata;
857 /* Tell our PME-only rank to switch grid */
858 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
863 print_grid(nullptr, debug, "", "switched to", set, -1);
866 if (pme_lb->stage == pme_lb->nstage)
868 print_grid(fp_err, fp_log, "", "optimal", set, -1);
872 /*! \brief Prepare for another round of PME load balancing
874 * \param[in,out] pme_lb Pointer to PME load balancing struct
875 * \param[in] bDlbUnlocked TRUE is DLB was locked and is now unlocked
877 * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
878 * the PP/PME balance might change and re-balancing can improve performance.
879 * This function adds 2 stages and adjusts the considered setup range.
881 static void continue_pme_loadbal(pme_load_balancing_t* pme_lb, gmx_bool bDlbUnlocked)
883 /* Add 2 tuning stages, keep the detected end of the setup range */
885 if (bDlbUnlocked && pme_lb->bSepPMERanks)
887 /* With separate PME ranks, DLB should always lower the PP load and
888 * can only increase the PME load (more communication and imbalance),
889 * so we only need to scan longer cut-off's.
891 pme_lb->lower_limit = pme_lb->cur;
893 pme_lb->start = pme_lb->lower_limit;
896 void pme_loadbal_do(pme_load_balancing_t* pme_lb,
900 const gmx::MDLogger& mdlog,
901 const t_inputrec& ir,
904 gmx::ArrayRef<const gmx::RVec> x,
905 gmx_wallcycle_t wcycle,
909 bool useGpuPmePpCommunication)
914 assert(pme_lb != nullptr);
916 if (!pme_lb->bActive)
921 n_prev = pme_lb->cycles_n;
922 cycles_prev = pme_lb->cycles_c;
923 wallcycle_get(wcycle, ewcSTEP, &pme_lb->cycles_n, &pme_lb->cycles_c);
925 /* Before the first step we haven't done any steps yet.
926 * Also handle cases where ir.init_step % ir.nstlist != 0.
927 * We also want to skip a number of steps and seconds while
928 * the CPU and GPU, when used, performance stabilizes.
930 if (!PAR(cr) || (DOMAINDECOMP(cr) && DDMASTER(cr->dd)))
932 pme_lb->startupTimeDelayElapsed = (gmx_gettime() - pme_lb->startTime < c_startupTimeDelay);
934 if (DOMAINDECOMP(cr))
936 dd_bcast(cr->dd, sizeof(bool), &pme_lb->startupTimeDelayElapsed);
939 if (pme_lb->cycles_n == 0 || step_rel < c_numFirstTuningIntervalSkip * ir.nstlist
940 || pme_lb->startupTimeDelayElapsed)
945 /* Sanity check, we expect nstlist cycle counts */
946 if (pme_lb->cycles_n - n_prev != ir.nstlist)
948 /* We could return here, but it's safer to issue an error and quit */
949 gmx_incons("pme_loadbal_do called at an interval != nstlist");
952 /* PME grid + cut-off optimization with GPUs or PME ranks */
953 if (!pme_lb->bBalance && pme_lb->bSepPMERanks)
955 if (pme_lb->bTriggerOnDLB)
957 pme_lb->bBalance = dd_dlb_is_on(cr->dd);
959 /* We should ignore the first timing to avoid timing allocation
960 * overhead. And since the PME load balancing is called just
961 * before DD repartitioning, the ratio returned by dd_pme_f_ratio
962 * is not over the last nstlist steps, but the nstlist steps before
963 * that. So the first useful ratio is available at step_rel=3*nstlist.
965 else if (step_rel >= c_numFirstTuningIntervalSkipWithSepPme * ir.nstlist)
967 GMX_ASSERT(DOMAINDECOMP(cr), "Domain decomposition should be active here");
968 if (DDMASTER(cr->dd))
970 /* If PME rank load is too high, start tuning. If
971 PME-PP direct GPU communication is active,
972 unconditionally start tuning since ratio will be
973 unreliable due to CPU-GPU asynchronicity in codepath */
974 pme_lb->bBalance = useGpuPmePpCommunication
976 : (dd_pme_f_ratio(cr->dd) >= loadBalanceTriggerFactor);
978 dd_bcast(cr->dd, sizeof(gmx_bool), &pme_lb->bBalance);
981 pme_lb->bActive = (pme_lb->bBalance || step_rel <= pme_lb->step_rel_stop);
984 /* The location in the code of this balancing termination is strange.
985 * You would expect to have it after the call to pme_load_balance()
986 * below, since there pme_lb->stage is updated.
987 * But when terminating directly after deciding on and selecting the
988 * optimal setup, DLB will turn on right away if it was locked before.
989 * This might be due to PME reinitialization. So we check stage here
990 * to allow for another nstlist steps with DLB locked to stabilize
993 if (pme_lb->bBalance && pme_lb->stage == pme_lb->nstage)
995 pme_lb->bBalance = FALSE;
997 if (DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
999 /* Unlock the DLB=auto, DLB is allowed to activate */
1000 dd_dlb_unlock(cr->dd);
1001 GMX_LOG(mdlog.warning)
1003 .appendText("NOTE: DLB can now turn on, when beneficial");
1005 /* We don't deactivate the tuning yet, since we will balance again
1006 * after DLB gets turned on, if it does within PMETune_period.
1008 continue_pme_loadbal(pme_lb, TRUE);
1009 pme_lb->bTriggerOnDLB = TRUE;
1010 pme_lb->step_rel_stop = step_rel + PMETunePeriod * ir.nstlist;
1014 /* We're completely done with PME tuning */
1015 pme_lb->bActive = FALSE;
1018 if (DOMAINDECOMP(cr))
1020 /* Set the cut-off limit to the final selected cut-off,
1021 * so we don't have artificial DLB limits.
1022 * This also ensures that we won't disable the currently
1023 * optimal setting during a second round of PME balancing.
1025 set_dd_dlb_max_cutoff(cr, fr->nbv->pairlistOuterRadius());
1029 if (pme_lb->bBalance)
1031 /* We might not have collected nstlist steps in cycles yet,
1032 * since init_step might not be a multiple of nstlist,
1033 * but the first data collected is skipped anyhow.
1035 pme_load_balance(pme_lb, cr, fp_err, fp_log, mdlog, ir, box, x,
1036 pme_lb->cycles_c - cycles_prev, fr->ic, fr->nbv.get(), &fr->pmedata, step);
1038 /* Update deprecated rlist in forcerec to stay in sync with fr->nbv */
1039 fr->rlist = fr->nbv->pairlistOuterRadius();
1041 if (ir.eDispCorr != edispcNO)
1043 fr->dispersionCorrection->setParameters(*fr->ic);
1047 if (!pme_lb->bBalance && (!pme_lb->bSepPMERanks || step_rel > pme_lb->step_rel_stop))
1049 /* We have just deactivated the balancing and we're not measuring PP/PME
1050 * imbalance during the first steps of the run: deactivate the tuning.
1052 pme_lb->bActive = FALSE;
1055 if (!(pme_lb->bActive) && DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
1057 /* Make sure DLB is allowed when we deactivate PME tuning */
1058 dd_dlb_unlock(cr->dd);
1059 GMX_LOG(mdlog.warning)
1061 .appendText("NOTE: DLB can now turn on, when beneficial");
1064 *bPrinting = pme_lb->bBalance;
1067 /*! \brief Return product of the number of PME grid points in each dimension */
1068 static int pme_grid_points(const pme_setup_t* setup)
1070 return setup->grid[XX] * setup->grid[YY] * setup->grid[ZZ];
1073 /*! \brief Print one load-balancing setting */
1074 static void print_pme_loadbal_setting(FILE* fplog, const char* name, const pme_setup_t* setup)
1076 fprintf(fplog, " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n", name,
1077 setup->rcut_coulomb, setup->rlistInner, setup->grid[XX], setup->grid[YY],
1078 setup->grid[ZZ], setup->spacing, 1 / setup->ewaldcoeff_q);
1081 /*! \brief Print all load-balancing settings */
1082 static void print_pme_loadbal_settings(pme_load_balancing_t* pme_lb,
1084 const gmx::MDLogger& mdlog,
1085 gmx_bool bNonBondedOnGPU)
1087 double pp_ratio, grid_ratio;
1088 real pp_ratio_temporary;
1090 pp_ratio_temporary = pme_lb->setup[pme_lb->cur].rlistInner / pme_lb->setup[0].rlistInner;
1091 pp_ratio = gmx::power3(pp_ratio_temporary);
1092 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])
1093 / static_cast<double>(pme_grid_points(&pme_lb->setup[0]));
1095 fprintf(fplog, "\n");
1096 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
1097 fprintf(fplog, "\n");
1098 /* Here we only warn when the optimal setting is the last one */
1099 if (pme_lb->elimited != epmelblimNO && pme_lb->cur == pme_loadbal_end(pme_lb) - 1)
1101 fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
1102 pmelblim_str[pme_lb->elimited]);
1103 fprintf(fplog, " you might not have reached a good load balance.\n");
1104 if (pme_lb->elimited == epmelblimDD)
1106 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
1108 fprintf(fplog, "\n");
1110 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
1111 fprintf(fplog, " particle-particle PME\n");
1112 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
1113 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
1114 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
1115 fprintf(fplog, " cost-ratio %4.2f %4.2f\n", pp_ratio, grid_ratio);
1116 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
1118 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
1120 GMX_LOG(mdlog.warning)
1123 "NOTE: PME load balancing increased the non-bonded workload by more than "
1125 " For better performance, use (more) PME ranks (mdrun -npme),\n"
1126 " or if you are beyond the scaling limit, use fewer total ranks (or "
1131 fprintf(fplog, "\n");
1135 void pme_loadbal_done(pme_load_balancing_t* pme_lb, FILE* fplog, const gmx::MDLogger& mdlog, gmx_bool bNonBondedOnGPU)
1137 if (fplog != nullptr && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
1139 print_pme_loadbal_settings(pme_lb, fplog, mdlog, bNonBondedOnGPU);