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, 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"
85 /*! \brief Parameters and settings for one PP-PME setup */
88 real rcut_coulomb; /**< Coulomb cut-off */
89 real rlistOuter; /**< cut-off for the outer pair-list */
90 real rlistInner; /**< cut-off for the inner pair-list */
91 real spacing; /**< (largest) PME grid spacing */
92 ivec grid; /**< the PME grid dimensions */
93 real grid_efficiency; /**< ineffiency factor for non-uniform grids <= 1 */
94 real ewaldcoeff_q; /**< Electrostatic Ewald coefficient */
95 real ewaldcoeff_lj; /**< LJ Ewald coefficient, only for the call to send_switchgrid */
96 struct gmx_pme_t* pmedata; /**< the data structure used in the PME code */
97 int count; /**< number of times this setup has been timed */
98 double cycles; /**< the fastest time for this setup in cycles */
101 /*! \brief After 50 nstlist periods of not observing imbalance: never tune PME */
102 const int PMETunePeriod = 50;
103 /*! \brief Trigger PME load balancing at more than 5% PME overload */
104 const real loadBalanceTriggerFactor = 1.05;
105 /*! \brief Scale the grid by a most at factor 1.7.
107 * This still leaves room for about 4-4.5x decrease in grid spacing while limiting the cases where
108 * large imbalance leads to extreme cutoff scaling for marginal benefits.
110 * This should help to avoid:
111 * - large increase in power consumption for little performance gain
112 * - increasing communication volume
115 const real c_maxSpacingScaling = 1.7;
116 /*! \brief In the initial scan, step by grids that are at least a factor 0.8 coarser */
117 const real gridpointsScaleFactor = 0.8;
118 /*! \brief In the initial scan, try to skip grids with uneven x/y/z spacing,
119 * checking if the "efficiency" is more than 5% worse than the previous grid.
121 const real relativeEfficiencyFactor = 1.05;
122 /*! \brief Rerun until a run is 12% slower setups than the fastest run so far */
123 const real maxRelativeSlowdownAccepted = 1.12;
124 /*! \brief If setups get more than 2% faster, do another round to avoid
125 * choosing a slower setup due to acceleration or fluctuations.
127 const real maxFluctuationAccepted = 1.02;
129 //! \brief Number of nstlist long tuning intervals to skip before starting
130 // load-balancing at the beginning of the run.
131 const int c_numFirstTuningIntervalSkip = 5;
132 //! \brief Number of nstlist long tuning intervals to skip before starting
133 // load-balancing at the beginning of the run with separate PME ranks. */
134 const int c_numFirstTuningIntervalSkipWithSepPme = 3;
135 //! \brief Number of nstlist long tuning intervals to skip after switching to a new setting
137 const int c_numPostSwitchTuningIntervalSkip = 1;
138 //! \brief Number of seconds to delay the tuning at startup to allow processors clocks to ramp up.
139 const double c_startupTimeDelay = 5.0;
141 /*! \brief Enumeration whose values describe the effect limiting the load balancing */
152 /*! \brief Descriptive strings matching ::epmelb */
153 static const char* pmelblim_str[epmelblimNR] = { "no", "box size", "domain decompostion",
154 "PME grid restriction",
155 "maximum allowed grid scaling" };
157 struct pme_load_balancing_t
159 gmx_bool bSepPMERanks; /**< do we have separate PME ranks? */
160 gmx_bool bActive; /**< is PME tuning active? */
161 int64_t step_rel_stop; /**< stop the tuning after this value of step_rel */
162 gmx_bool bTriggerOnDLB; /**< trigger balancing only on DD DLB */
163 gmx_bool bBalance; /**< are we in the balancing phase, i.e. trying different setups? */
164 int nstage; /**< the current maximum number of stages */
165 bool startupTimeDelayElapsed; /**< Has the c_startupTimeDelay elapsed indicating that the balancing can start. */
167 real cut_spacing; /**< the minimum cutoff / PME grid spacing ratio */
168 real rcut_vdw; /**< Vdw cutoff (does not change) */
169 real rcut_coulomb_start; /**< Initial electrostatics cutoff */
170 real rbufOuter_coulomb; /**< the outer pairlist buffer size */
171 real rbufOuter_vdw; /**< the outer pairlist buffer size */
172 real rbufInner_coulomb; /**< the inner pairlist buffer size */
173 real rbufInner_vdw; /**< the inner pairlist buffer size */
174 matrix box_start; /**< the initial simulation box */
175 std::vector<pme_setup_t> setup; /**< the PME+cutoff setups */
176 int cur; /**< the index (in setup) of the current setup */
177 int fastest; /**< index of the fastest setup up till now */
178 int lower_limit; /**< don't go below this setup index */
179 int start; /**< start of setup index range to consider in stage>0 */
180 int end; /**< end of setup index range to consider in stage>0 */
181 int elimited; /**< was the balancing limited, uses enum above */
182 int cutoff_scheme; /**< Verlet or group cut-offs */
184 int stage; /**< the current stage */
186 int cycles_n; /**< step cycle counter cumulative count */
187 double cycles_c; /**< step cycle counter cumulative cycles */
188 double startTime; /**< time stamp when the balancing was started on the master rank (relative to the UNIX epoch start).*/
191 /* TODO The code in this file should call this getter, rather than
192 * read bActive anywhere */
193 bool pme_loadbal_is_active(const pme_load_balancing_t* pme_lb)
195 return pme_lb != nullptr && pme_lb->bActive;
198 // TODO Return a unique_ptr to pme_load_balancing_t
199 void pme_loadbal_init(pme_load_balancing_t** pme_lb_p,
201 const gmx::MDLogger& mdlog,
202 const t_inputrec& ir,
204 const interaction_const_t& ic,
205 const nonbonded_verlet_t& nbv,
210 pme_load_balancing_t* pme_lb;
214 // Note that we don't (yet) support PME load balancing with LJ-PME only.
215 GMX_RELEASE_ASSERT(EEL_PME(ir.coulombtype),
216 "pme_loadbal_init called without PME electrostatics");
217 // To avoid complexity, we require a single cut-off with PME for q+LJ.
218 // This is checked by grompp, but it doesn't hurt to check again.
219 GMX_RELEASE_ASSERT(!(EEL_PME(ir.coulombtype) && EVDW_PME(ir.vdwtype) && ir.rcoulomb != ir.rvdw),
220 "With Coulomb and LJ PME, rcoulomb should be equal to rvdw");
222 pme_lb = new pme_load_balancing_t;
224 pme_lb->bSepPMERanks = !thisRankHasDuty(cr, DUTY_PME);
226 /* Initially we turn on balancing directly on based on PP/PME imbalance */
227 pme_lb->bTriggerOnDLB = FALSE;
229 /* Any number of stages >= 2 is supported */
232 pme_lb->cutoff_scheme = ir.cutoff_scheme;
234 pme_lb->rbufOuter_coulomb = nbv.pairlistOuterRadius() - ic.rcoulomb;
235 pme_lb->rbufOuter_vdw = nbv.pairlistOuterRadius() - ic.rvdw;
236 pme_lb->rbufInner_coulomb = nbv.pairlistInnerRadius() - ic.rcoulomb;
237 pme_lb->rbufInner_vdw = nbv.pairlistInnerRadius() - ic.rvdw;
239 /* Scale box with Ewald wall factor; note that we pmedata->boxScaler
240 * can't always usedd as it's not available with separate PME ranks.
242 EwaldBoxZScaler boxScaler(ir);
243 boxScaler.scaleBox(box, pme_lb->box_start);
245 pme_lb->setup.resize(1);
247 pme_lb->rcut_vdw = ic.rvdw;
248 pme_lb->rcut_coulomb_start = ir.rcoulomb;
251 pme_lb->setup[0].rcut_coulomb = ic.rcoulomb;
252 pme_lb->setup[0].rlistOuter = nbv.pairlistOuterRadius();
253 pme_lb->setup[0].rlistInner = nbv.pairlistInnerRadius();
254 pme_lb->setup[0].grid[XX] = ir.nkx;
255 pme_lb->setup[0].grid[YY] = ir.nky;
256 pme_lb->setup[0].grid[ZZ] = ir.nkz;
257 pme_lb->setup[0].ewaldcoeff_q = ic.ewaldcoeff_q;
258 pme_lb->setup[0].ewaldcoeff_lj = ic.ewaldcoeff_lj;
260 if (!pme_lb->bSepPMERanks)
262 GMX_RELEASE_ASSERT(pmedata,
263 "On ranks doing both PP and PME we need a valid pmedata object");
264 pme_lb->setup[0].pmedata = pmedata;
268 for (d = 0; d < DIM; d++)
270 sp = norm(pme_lb->box_start[d]) / pme_lb->setup[0].grid[d];
276 pme_lb->setup[0].spacing = spm;
278 if (ir.fourier_spacing > 0)
280 pme_lb->cut_spacing = ir.rcoulomb / ir.fourier_spacing;
284 pme_lb->cut_spacing = ir.rcoulomb / pme_lb->setup[0].spacing;
290 pme_lb->lower_limit = 0;
293 pme_lb->elimited = epmelblimNO;
295 pme_lb->cycles_n = 0;
296 pme_lb->cycles_c = 0;
297 // only master ranks do timing
298 if (!PAR(cr) || (DOMAINDECOMP(cr) && DDMASTER(cr->dd)))
300 pme_lb->startTime = gmx_gettime();
303 if (!wallcycle_have_counter())
305 GMX_LOG(mdlog.warning)
308 "NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use "
309 "PME-PP balancing.");
312 /* Tune with GPUs and/or separate PME ranks.
313 * When running only on a CPU without PME ranks, PME tuning will only help
314 * with small numbers of atoms in the cut-off sphere.
316 pme_lb->bActive = (wallcycle_have_counter() && (bUseGPU || pme_lb->bSepPMERanks));
318 /* With GPUs and no separate PME ranks we can't measure the PP/PME
319 * imbalance, so we start balancing right away.
320 * Otherwise we only start balancing after we observe imbalance.
322 pme_lb->bBalance = (pme_lb->bActive && (bUseGPU && !pme_lb->bSepPMERanks));
324 pme_lb->step_rel_stop = PMETunePeriod * ir.nstlist;
326 /* Delay DD load balancing when GPUs are used */
327 if (pme_lb->bActive && DOMAINDECOMP(cr) && cr->dd->nnodes > 1 && bUseGPU)
329 /* Lock DLB=auto to off (does nothing when DLB=yes/no.
330 * With GPUs and separate PME nodes, we want to first
331 * do PME tuning without DLB, since DLB might limit
332 * the cut-off, which never improves performance.
333 * We allow for DLB + PME tuning after a first round of tuning.
336 if (dd_dlb_is_locked(cr->dd))
338 GMX_LOG(mdlog.warning)
340 .appendText("NOTE: DLB will not turn on during the first phase of PME tuning");
347 /*! \brief Try to increase the cutoff during load balancing */
348 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t* pme_lb, int pme_order, const gmx_domdec_t* dd)
351 real tmpr_coulomb, tmpr_vdw;
355 /* Try to add a new setup with next larger cut-off to the list */
358 set.pmedata = nullptr;
360 NumPmeDomains numPmeDomains = getNumPmeDomains(dd);
365 /* Avoid infinite while loop, which can occur at the minimum grid size.
366 * Note that in practice load balancing will stop before this point.
367 * The factor 2.1 allows for the extreme case in which only grids
368 * of powers of 2 are allowed (the current code supports more grids).
376 clear_ivec(set.grid);
377 sp = calcFftGrid(nullptr, pme_lb->box_start, fac * pme_lb->setup[pme_lb->cur].spacing,
378 minimalPmeGridSize(pme_order), &set.grid[XX], &set.grid[YY], &set.grid[ZZ]);
380 /* As here we can't easily check if one of the PME ranks
381 * uses threading, we do a conservative grid check.
382 * This means we can't use pme_order or less grid lines
383 * per PME rank along x, which is not a strong restriction.
385 grid_ok = gmx_pme_check_restrictions(pme_order, set.grid[XX], set.grid[YY], set.grid[ZZ],
386 numPmeDomains.x, true, false);
387 } while (sp <= 1.001 * pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
389 set.rcut_coulomb = pme_lb->cut_spacing * sp;
390 if (set.rcut_coulomb < pme_lb->rcut_coulomb_start)
392 /* This is unlikely, but can happen when e.g. continuing from
393 * a checkpoint after equilibration where the box shrank a lot.
394 * We want to avoid rcoulomb getting smaller than rvdw
395 * and there might be more issues with decreasing rcoulomb.
397 set.rcut_coulomb = pme_lb->rcut_coulomb_start;
400 if (pme_lb->cutoff_scheme == ecutsVERLET)
402 /* Never decrease the Coulomb and VdW list buffers */
403 set.rlistOuter = std::max(set.rcut_coulomb + pme_lb->rbufOuter_coulomb,
404 pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw);
405 set.rlistInner = std::max(set.rcut_coulomb + pme_lb->rbufInner_coulomb,
406 pme_lb->rcut_vdw + pme_lb->rbufInner_vdw);
410 /* TODO Remove these lines and pme_lb->cutoff_scheme */
411 tmpr_coulomb = set.rcut_coulomb + pme_lb->rbufOuter_coulomb;
412 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw;
413 /* Two (known) bugs with cutoff-scheme=group here:
414 * - This modification of rlist results in incorrect DD comunication.
415 * - We should set fr->bTwinRange = (fr->rlistlong > fr->rlist).
417 set.rlistOuter = std::min(tmpr_coulomb, tmpr_vdw);
418 set.rlistInner = set.rlistOuter;
422 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
423 set.grid_efficiency = 1;
424 for (d = 0; d < DIM; d++)
426 set.grid_efficiency *= (set.grid[d] * sp) / norm(pme_lb->box_start[d]);
428 /* The Ewald coefficient is inversly proportional to the cut-off */
429 set.ewaldcoeff_q = pme_lb->setup[0].ewaldcoeff_q * pme_lb->setup[0].rcut_coulomb / set.rcut_coulomb;
430 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
431 set.ewaldcoeff_lj = pme_lb->setup[0].ewaldcoeff_lj * pme_lb->setup[0].rcut_coulomb / set.rcut_coulomb;
438 fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n", set.grid[XX],
439 set.grid[YY], set.grid[ZZ], set.rcut_coulomb);
441 pme_lb->setup.push_back(set);
445 /*! \brief Print the PME grid */
446 static void print_grid(FILE* fp_err, FILE* fp_log, const char* pre, const char* desc, const pme_setup_t* set, double cycles)
448 auto buf = gmx::formatString("%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f", pre, desc,
449 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
452 buf += gmx::formatString(": %.1f M-cycles", cycles * 1e-6);
454 if (fp_err != nullptr)
456 fprintf(fp_err, "\r%s\n", buf.c_str());
459 if (fp_log != nullptr)
461 fprintf(fp_log, "%s\n", buf.c_str());
465 /*! \brief Return the index of the last setup used in PME load balancing */
466 static int pme_loadbal_end(pme_load_balancing_t* pme_lb)
468 /* In the initial stage only n is set; end is not set yet */
475 return pme_lb->setup.size();
479 /*! \brief Print descriptive string about what limits PME load balancing */
480 static void print_loadbal_limited(FILE* fp_err, FILE* fp_log, int64_t step, pme_load_balancing_t* pme_lb)
482 auto buf = gmx::formatString(
483 "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
484 gmx::int64ToString(step).c_str(), pmelblim_str[pme_lb->elimited],
485 pme_lb->setup[pme_loadbal_end(pme_lb) - 1].rcut_coulomb);
486 if (fp_err != nullptr)
488 fprintf(fp_err, "\r%s\n", buf.c_str());
491 if (fp_log != nullptr)
493 fprintf(fp_log, "%s\n", buf.c_str());
497 /*! \brief Switch load balancing to stage 1
499 * In this stage, only reasonably fast setups are run again. */
500 static void switch_to_stage1(pme_load_balancing_t* pme_lb)
502 /* Increase start until we find a setup that is not slower than
503 * maxRelativeSlowdownAccepted times the fastest setup.
505 pme_lb->start = pme_lb->lower_limit;
506 while (pme_lb->start + 1 < gmx::ssize(pme_lb->setup)
507 && (pme_lb->setup[pme_lb->start].count == 0
508 || pme_lb->setup[pme_lb->start].cycles
509 > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted))
513 /* While increasing start, we might have skipped setups that we did not
514 * time during stage 0. We want to extend the range for stage 1 to include
515 * any skipped setups that lie between setups that were measured to be
516 * acceptably fast and too slow.
518 while (pme_lb->start > pme_lb->lower_limit && pme_lb->setup[pme_lb->start - 1].count == 0)
523 /* Decrease end only with setups that we timed and that are slow. */
524 pme_lb->end = pme_lb->setup.size();
525 if (pme_lb->setup[pme_lb->end - 1].count > 0
526 && pme_lb->setup[pme_lb->end - 1].cycles
527 > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted)
534 /* Next we want to choose setup pme_lb->end-1, but as we will decrease
535 * pme_lb->cur by one right after returning, we set cur to end.
537 pme_lb->cur = pme_lb->end;
540 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
542 * The adjustment is done to generate a different non-bonded PP and PME load.
543 * With separate PME ranks (PP and PME on different processes) or with
544 * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
545 * and changing the load will affect the load balance and performance.
546 * The total time for a set of integration steps is monitored and a range
547 * of grid/cut-off setups is scanned. After calling pme_load_balance many
548 * times and acquiring enough statistics, the best performing setup is chosen.
549 * Here we try to take into account fluctuations and changes due to external
550 * factors as well as DD load balancing.
552 static void pme_load_balance(pme_load_balancing_t* pme_lb,
556 const gmx::MDLogger& mdlog,
557 const t_inputrec& ir,
559 gmx::ArrayRef<const gmx::RVec> x,
561 interaction_const_t* ic,
562 struct nonbonded_verlet_t* nbv,
563 struct gmx_pme_t** pmedata,
569 char buf[STRLEN], sbuf[22];
573 gmx_sumd(1, &cycles, cr);
574 cycles /= cr->nnodes;
577 set = &pme_lb->setup[pme_lb->cur];
580 /* Skip the first c_numPostSwitchTuningIntervalSkip cycles because the first step
581 * after a switch is much slower due to allocation and/or caching effects.
583 if (set->count % (c_numPostSwitchTuningIntervalSkip + 1) != 0)
588 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
589 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
591 GMX_RELEASE_ASSERT(set->count > c_numPostSwitchTuningIntervalSkip, "We should skip cycles");
592 if (set->count == (c_numPostSwitchTuningIntervalSkip + 1))
594 set->cycles = cycles;
598 if (cycles * maxFluctuationAccepted < set->cycles && pme_lb->stage == pme_lb->nstage - 1)
600 /* The performance went up a lot (due to e.g. DD load balancing).
601 * Add a stage, keep the minima, but rescan all setups.
608 "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this "
610 "Increased the number stages to %d"
611 " and ignoring the previous performance\n",
612 set->grid[XX], set->grid[YY], set->grid[ZZ], set->cycles * 1e-6,
613 cycles * 1e-6, maxFluctuationAccepted, pme_lb->nstage);
616 set->cycles = std::min(set->cycles, cycles);
619 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
621 pme_lb->fastest = pme_lb->cur;
623 if (DOMAINDECOMP(cr))
625 /* We found a new fastest setting, ensure that with subsequent
626 * shorter cut-off's the dynamic load balancing does not make
627 * the use of the current cut-off impossible. This solution is
628 * a trade-off, as the PME load balancing and DD domain size
629 * load balancing can interact in complex ways.
630 * With the Verlet kernels, DD load imbalance will usually be
631 * mainly due to bonded interaction imbalance, which will often
632 * quickly push the domain boundaries beyond the limit for the
633 * optimal, PME load balanced, cut-off. But it could be that
634 * better overal performance can be obtained with a slightly
635 * shorter cut-off and better DD load balancing.
637 set_dd_dlb_max_cutoff(cr, pme_lb->setup[pme_lb->fastest].rlistOuter);
640 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
642 /* Check in stage 0 if we should stop scanning grids.
643 * Stop when the time is more than maxRelativeSlowDownAccepted longer than the fastest.
645 if (pme_lb->stage == 0 && pme_lb->cur > 0
646 && cycles > pme_lb->setup[pme_lb->fastest].cycles * maxRelativeSlowdownAccepted)
648 pme_lb->setup.resize(pme_lb->cur + 1);
649 /* Done with scanning, go to stage 1 */
650 switch_to_stage1(pme_lb);
653 if (pme_lb->stage == 0)
657 gridsize_start = set->grid[XX] * set->grid[YY] * set->grid[ZZ];
661 if (pme_lb->cur + 1 < gmx::ssize(pme_lb->setup))
663 /* We had already generated the next setup */
668 /* Find the next setup */
669 OK = pme_loadbal_increase_cutoff(pme_lb, ir.pme_order, cr->dd);
673 pme_lb->elimited = epmelblimPMEGRID;
678 && pme_lb->setup[pme_lb->cur + 1].spacing > c_maxSpacingScaling * pme_lb->setup[0].spacing)
681 pme_lb->elimited = epmelblimMAXSCALING;
684 if (OK && ir.pbcType != PbcType::No)
686 OK = (gmx::square(pme_lb->setup[pme_lb->cur + 1].rlistOuter)
687 <= max_cutoff2(ir.pbcType, box));
690 pme_lb->elimited = epmelblimBOX;
698 if (DOMAINDECOMP(cr))
700 OK = change_dd_cutoff(cr, box, x, pme_lb->setup[pme_lb->cur].rlistOuter);
703 /* Failed: do not use this setup */
705 pme_lb->elimited = epmelblimDD;
711 /* We hit the upper limit for the cut-off,
712 * the setup should not go further than cur.
714 pme_lb->setup.resize(pme_lb->cur + 1);
715 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
716 /* Switch to the next stage */
717 switch_to_stage1(pme_lb);
720 && !(pme_lb->setup[pme_lb->cur].grid[XX] * pme_lb->setup[pme_lb->cur].grid[YY]
721 * pme_lb->setup[pme_lb->cur].grid[ZZ]
722 < gridsize_start * gridpointsScaleFactor
723 && pme_lb->setup[pme_lb->cur].grid_efficiency
724 < pme_lb->setup[pme_lb->cur - 1].grid_efficiency * relativeEfficiencyFactor));
727 if (pme_lb->stage > 0 && pme_lb->end == 1)
729 pme_lb->cur = pme_lb->lower_limit;
730 pme_lb->stage = pme_lb->nstage;
732 else if (pme_lb->stage > 0 && pme_lb->end > 1)
734 /* If stage = nstage-1:
735 * scan over all setups, rerunning only those setups
736 * which are not much slower than the fastest
739 * Note that we loop backward to minimize the risk of the cut-off
740 * getting limited by DD DLB, since the DLB cut-off limit is set
741 * to the fastest PME setup.
745 if (pme_lb->cur > pme_lb->start)
753 pme_lb->cur = pme_lb->end - 1;
755 } while (pme_lb->stage == pme_lb->nstage - 1 && pme_lb->setup[pme_lb->cur].count > 0
756 && pme_lb->setup[pme_lb->cur].cycles > cycles_fast * maxRelativeSlowdownAccepted);
758 if (pme_lb->stage == pme_lb->nstage)
760 /* We are done optimizing, use the fastest setup we found */
761 pme_lb->cur = pme_lb->fastest;
765 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
767 OK = change_dd_cutoff(cr, box, x, pme_lb->setup[pme_lb->cur].rlistOuter);
770 /* For some reason the chosen cut-off is incompatible with DD.
771 * We should continue scanning a more limited range of cut-off's.
773 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
775 /* stage=nstage says we're finished, but we should continue
776 * balancing, so we set back stage which was just incremented.
780 if (pme_lb->cur <= pme_lb->fastest)
782 /* This should not happen, as we set limits on the DLB bounds.
783 * But we implement a complete failsafe solution anyhow.
785 GMX_LOG(mdlog.warning)
787 .appendTextFormatted(
788 "The fastest PP/PME load balancing setting (cutoff %.3d nm) is no "
789 "longer available due to DD DLB or box size limitations",
791 pme_lb->fastest = pme_lb->lower_limit;
792 pme_lb->start = pme_lb->lower_limit;
794 /* Limit the range to below the current cut-off, scan from start */
795 pme_lb->end = pme_lb->cur;
796 pme_lb->cur = pme_lb->start;
797 pme_lb->elimited = epmelblimDD;
798 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
802 /* Change the Coulomb cut-off and the PME grid */
804 set = &pme_lb->setup[pme_lb->cur];
806 ic->rcoulomb = set->rcut_coulomb;
807 nbv->changePairlistRadii(set->rlistOuter, set->rlistInner);
808 ic->ewaldcoeff_q = set->ewaldcoeff_q;
809 /* TODO: centralize the code that sets the potentials shifts */
810 if (ic->coulomb_modifier == eintmodPOTSHIFT)
812 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
813 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
815 if (EVDW_PME(ic->vdwtype))
817 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
818 ic->rvdw = set->rcut_coulomb;
819 ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
820 if (ic->vdw_modifier == eintmodPOTSHIFT)
824 ic->dispersion_shift.cpot = -1.0 / gmx::power6(static_cast<double>(ic->rvdw));
825 ic->repulsion_shift.cpot = -1.0 / gmx::power12(static_cast<double>(ic->rvdw));
826 crc2 = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
828 (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
832 /* We always re-initialize the tables whether they are used or not */
833 init_interaction_const_tables(nullptr, ic);
835 Nbnxm::gpu_pme_loadbal_update_param(nbv, ic);
837 if (!pme_lb->bSepPMERanks)
840 * CPU PME keeps a list of allocated pmedata's, that's why pme_lb->setup[pme_lb->cur].pmedata is not always nullptr.
841 * GPU PME, however, currently needs the gmx_pme_reinit always called on load balancing
842 * (pme_gpu_reinit might be not sufficiently decoupled from gmx_pme_init).
843 * This can lead to a lot of reallocations for PME GPU.
844 * Would be nicer if the allocated grid list was hidden within a single pmedata structure.
846 if ((pme_lb->setup[pme_lb->cur].pmedata == nullptr)
847 || pme_gpu_task_enabled(pme_lb->setup[pme_lb->cur].pmedata))
849 /* Generate a new PME data structure,
850 * copying part of the old pointers.
852 gmx_pme_reinit(&set->pmedata, cr, pme_lb->setup[0].pmedata, &ir, set->grid,
853 set->ewaldcoeff_q, set->ewaldcoeff_lj);
855 *pmedata = set->pmedata;
859 /* Tell our PME-only rank to switch grid */
860 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
865 print_grid(nullptr, debug, "", "switched to", set, -1);
868 if (pme_lb->stage == pme_lb->nstage)
870 print_grid(fp_err, fp_log, "", "optimal", set, -1);
874 /*! \brief Prepare for another round of PME load balancing
876 * \param[in,out] pme_lb Pointer to PME load balancing struct
877 * \param[in] bDlbUnlocked TRUE is DLB was locked and is now unlocked
879 * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
880 * the PP/PME balance might change and re-balancing can improve performance.
881 * This function adds 2 stages and adjusts the considered setup range.
883 static void continue_pme_loadbal(pme_load_balancing_t* pme_lb, gmx_bool bDlbUnlocked)
885 /* Add 2 tuning stages, keep the detected end of the setup range */
887 if (bDlbUnlocked && pme_lb->bSepPMERanks)
889 /* With separate PME ranks, DLB should always lower the PP load and
890 * can only increase the PME load (more communication and imbalance),
891 * so we only need to scan longer cut-off's.
893 pme_lb->lower_limit = pme_lb->cur;
895 pme_lb->start = pme_lb->lower_limit;
898 void pme_loadbal_do(pme_load_balancing_t* pme_lb,
902 const gmx::MDLogger& mdlog,
903 const t_inputrec& ir,
906 gmx::ArrayRef<const gmx::RVec> x,
907 gmx_wallcycle_t wcycle,
911 bool useGpuPmePpCommunication)
916 assert(pme_lb != nullptr);
918 if (!pme_lb->bActive)
923 n_prev = pme_lb->cycles_n;
924 cycles_prev = pme_lb->cycles_c;
925 wallcycle_get(wcycle, ewcSTEP, &pme_lb->cycles_n, &pme_lb->cycles_c);
927 /* Before the first step we haven't done any steps yet.
928 * Also handle cases where ir.init_step % ir.nstlist != 0.
929 * We also want to skip a number of steps and seconds while
930 * the CPU and GPU, when used, performance stabilizes.
932 if (!PAR(cr) || (DOMAINDECOMP(cr) && DDMASTER(cr->dd)))
934 pme_lb->startupTimeDelayElapsed = (gmx_gettime() - pme_lb->startTime < c_startupTimeDelay);
936 if (DOMAINDECOMP(cr))
938 dd_bcast(cr->dd, sizeof(bool), &pme_lb->startupTimeDelayElapsed);
941 if (pme_lb->cycles_n == 0 || step_rel < c_numFirstTuningIntervalSkip * ir.nstlist
942 || pme_lb->startupTimeDelayElapsed)
947 /* Sanity check, we expect nstlist cycle counts */
948 if (pme_lb->cycles_n - n_prev != ir.nstlist)
950 /* We could return here, but it's safer to issue an error and quit */
951 gmx_incons("pme_loadbal_do called at an interval != nstlist");
954 /* PME grid + cut-off optimization with GPUs or PME ranks */
955 if (!pme_lb->bBalance && pme_lb->bSepPMERanks)
957 if (pme_lb->bTriggerOnDLB)
959 pme_lb->bBalance = dd_dlb_is_on(cr->dd);
961 /* We should ignore the first timing to avoid timing allocation
962 * overhead. And since the PME load balancing is called just
963 * before DD repartitioning, the ratio returned by dd_pme_f_ratio
964 * is not over the last nstlist steps, but the nstlist steps before
965 * that. So the first useful ratio is available at step_rel=3*nstlist.
967 else if (step_rel >= c_numFirstTuningIntervalSkipWithSepPme * ir.nstlist)
969 GMX_ASSERT(DOMAINDECOMP(cr), "Domain decomposition should be active here");
970 if (DDMASTER(cr->dd))
972 /* If PME rank load is too high, start tuning. If
973 PME-PP direct GPU communication is active,
974 unconditionally start tuning since ratio will be
975 unreliable due to CPU-GPU asynchronicity in codepath */
976 pme_lb->bBalance = useGpuPmePpCommunication
978 : (dd_pme_f_ratio(cr->dd) >= loadBalanceTriggerFactor);
980 dd_bcast(cr->dd, sizeof(gmx_bool), &pme_lb->bBalance);
983 pme_lb->bActive = (pme_lb->bBalance || step_rel <= pme_lb->step_rel_stop);
986 /* The location in the code of this balancing termination is strange.
987 * You would expect to have it after the call to pme_load_balance()
988 * below, since there pme_lb->stage is updated.
989 * But when terminating directly after deciding on and selecting the
990 * optimal setup, DLB will turn on right away if it was locked before.
991 * This might be due to PME reinitialization. So we check stage here
992 * to allow for another nstlist steps with DLB locked to stabilize
995 if (pme_lb->bBalance && pme_lb->stage == pme_lb->nstage)
997 pme_lb->bBalance = FALSE;
999 if (DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
1001 /* Unlock the DLB=auto, DLB is allowed to activate */
1002 dd_dlb_unlock(cr->dd);
1003 GMX_LOG(mdlog.warning)
1005 .appendText("NOTE: DLB can now turn on, when beneficial");
1007 /* We don't deactivate the tuning yet, since we will balance again
1008 * after DLB gets turned on, if it does within PMETune_period.
1010 continue_pme_loadbal(pme_lb, TRUE);
1011 pme_lb->bTriggerOnDLB = TRUE;
1012 pme_lb->step_rel_stop = step_rel + PMETunePeriod * ir.nstlist;
1016 /* We're completely done with PME tuning */
1017 pme_lb->bActive = FALSE;
1020 if (DOMAINDECOMP(cr))
1022 /* Set the cut-off limit to the final selected cut-off,
1023 * so we don't have artificial DLB limits.
1024 * This also ensures that we won't disable the currently
1025 * optimal setting during a second round of PME balancing.
1027 set_dd_dlb_max_cutoff(cr, fr->nbv->pairlistOuterRadius());
1031 if (pme_lb->bBalance)
1033 /* We might not have collected nstlist steps in cycles yet,
1034 * since init_step might not be a multiple of nstlist,
1035 * but the first data collected is skipped anyhow.
1037 pme_load_balance(pme_lb, cr, fp_err, fp_log, mdlog, ir, box, x,
1038 pme_lb->cycles_c - cycles_prev, fr->ic, fr->nbv.get(), &fr->pmedata, step);
1040 /* Update deprecated rlist in forcerec to stay in sync with fr->nbv */
1041 fr->rlist = fr->nbv->pairlistOuterRadius();
1043 if (ir.eDispCorr != edispcNO)
1045 fr->dispersionCorrection->setParameters(*fr->ic);
1049 if (!pme_lb->bBalance && (!pme_lb->bSepPMERanks || step_rel > pme_lb->step_rel_stop))
1051 /* We have just deactivated the balancing and we're not measuring PP/PME
1052 * imbalance during the first steps of the run: deactivate the tuning.
1054 pme_lb->bActive = FALSE;
1057 if (!(pme_lb->bActive) && DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
1059 /* Make sure DLB is allowed when we deactivate PME tuning */
1060 dd_dlb_unlock(cr->dd);
1061 GMX_LOG(mdlog.warning)
1063 .appendText("NOTE: DLB can now turn on, when beneficial");
1066 *bPrinting = pme_lb->bBalance;
1069 /*! \brief Return product of the number of PME grid points in each dimension */
1070 static int pme_grid_points(const pme_setup_t* setup)
1072 return setup->grid[XX] * setup->grid[YY] * setup->grid[ZZ];
1075 /*! \brief Print one load-balancing setting */
1076 static void print_pme_loadbal_setting(FILE* fplog, const char* name, const pme_setup_t* setup)
1078 fprintf(fplog, " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n", name,
1079 setup->rcut_coulomb, setup->rlistInner, setup->grid[XX], setup->grid[YY],
1080 setup->grid[ZZ], setup->spacing, 1 / setup->ewaldcoeff_q);
1083 /*! \brief Print all load-balancing settings */
1084 static void print_pme_loadbal_settings(pme_load_balancing_t* pme_lb,
1086 const gmx::MDLogger& mdlog,
1087 gmx_bool bNonBondedOnGPU)
1089 double pp_ratio, grid_ratio;
1090 real pp_ratio_temporary;
1092 pp_ratio_temporary = pme_lb->setup[pme_lb->cur].rlistInner / pme_lb->setup[0].rlistInner;
1093 pp_ratio = gmx::power3(pp_ratio_temporary);
1094 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])
1095 / static_cast<double>(pme_grid_points(&pme_lb->setup[0]));
1097 fprintf(fplog, "\n");
1098 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
1099 fprintf(fplog, "\n");
1100 /* Here we only warn when the optimal setting is the last one */
1101 if (pme_lb->elimited != epmelblimNO && pme_lb->cur == pme_loadbal_end(pme_lb) - 1)
1103 fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
1104 pmelblim_str[pme_lb->elimited]);
1105 fprintf(fplog, " you might not have reached a good load balance.\n");
1106 if (pme_lb->elimited == epmelblimDD)
1108 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
1110 fprintf(fplog, "\n");
1112 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
1113 fprintf(fplog, " particle-particle PME\n");
1114 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
1115 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
1116 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
1117 fprintf(fplog, " cost-ratio %4.2f %4.2f\n", pp_ratio, grid_ratio);
1118 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
1120 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
1122 GMX_LOG(mdlog.warning)
1125 "NOTE: PME load balancing increased the non-bonded workload by more than "
1127 " For better performance, use (more) PME ranks (mdrun -npme),\n"
1128 " or if you are beyond the scaling limit, use fewer total ranks (or "
1133 fprintf(fplog, "\n");
1137 void pme_loadbal_done(pme_load_balancing_t* pme_lb, FILE* fplog, const gmx::MDLogger& mdlog, gmx_bool bNonBondedOnGPU)
1139 if (fplog != nullptr && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
1141 print_pme_loadbal_settings(pme_lb, fplog, mdlog, bNonBondedOnGPU);