2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * \brief This file contains function definitions necessary for
38 * managing automatic load balance of PME calculations (Coulomb and
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_ewald
46 #include "pme-load-balancing.h"
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_network.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/ewald/ewald-utils.h"
58 #include "gromacs/ewald/pme.h"
59 #include "gromacs/fft/calcgrid.h"
60 #include "gromacs/gmxlib/network.h"
61 #include "gromacs/math/functions.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/forcerec.h"
64 #include "gromacs/mdlib/nb_verlet.h"
65 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
66 #include "gromacs/mdlib/nbnxn_pairlist.h"
67 #include "gromacs/mdlib/sim_util.h"
68 #include "gromacs/mdtypes/commrec.h"
69 #include "gromacs/mdtypes/inputrec.h"
70 #include "gromacs/mdtypes/md_enums.h"
71 #include "gromacs/mdtypes/state.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/timing/wallcycle.h"
74 #include "gromacs/utility/cstringutil.h"
75 #include "gromacs/utility/fatalerror.h"
76 #include "gromacs/utility/gmxassert.h"
77 #include "gromacs/utility/logger.h"
78 #include "gromacs/utility/smalloc.h"
80 #include "pme-internal.h"
82 /*! \brief Parameters and settings for one PP-PME setup */
84 real rcut_coulomb; /**< Coulomb cut-off */
85 real rlistOuter; /**< cut-off for the outer pair-list */
86 real rlistInner; /**< cut-off for the inner pair-list */
87 real spacing; /**< (largest) PME grid spacing */
88 ivec grid; /**< the PME grid dimensions */
89 real grid_efficiency; /**< ineffiency factor for non-uniform grids <= 1 */
90 real ewaldcoeff_q; /**< Electrostatic Ewald coefficient */
91 real ewaldcoeff_lj; /**< LJ Ewald coefficient, only for the call to send_switchgrid */
92 struct gmx_pme_t *pmedata; /**< the data structure used in the PME code */
93 int count; /**< number of times this setup has been timed */
94 double cycles; /**< the fastest time for this setup in cycles */
97 /*! \brief After 50 nstlist periods of not observing imbalance: never tune PME */
98 const int PMETunePeriod = 50;
99 /*! \brief Trigger PME load balancing at more than 5% PME overload */
100 const real loadBalanceTriggerFactor = 1.05;
101 /*! \brief Scale the grid by a most at factor 1.7.
103 * This still leaves room for about 4-4.5x decrease in grid spacing while limiting the cases where
104 * large imbalance leads to extreme cutoff scaling for marginal benefits.
106 * This should help to avoid:
107 * - large increase in power consumption for little performance gain
108 * - increasing communication volume
111 const real c_maxSpacingScaling = 1.7;
112 /*! \brief In the initial scan, step by grids that are at least a factor 0.8 coarser */
113 const real gridpointsScaleFactor = 0.8;
114 /*! \brief In the initial scan, try to skip grids with uneven x/y/z spacing,
115 * checking if the "efficiency" is more than 5% worse than the previous grid.
117 const real relativeEfficiencyFactor = 1.05;
118 /*! \brief Rerun until a run is 12% slower setups than the fastest run so far */
119 const real maxRelativeSlowdownAccepted = 1.12;
120 /*! \brief If setups get more than 2% faster, do another round to avoid
121 * choosing a slower setup due to acceleration or fluctuations.
123 const real maxFluctuationAccepted = 1.02;
125 /*! \brief Enumeration whose values describe the effect limiting the load balancing */
127 epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimPMEGRID, epmelblimMAXSCALING, epmelblimNR
130 /*! \brief Descriptive strings matching ::epmelb */
131 static const char *pmelblim_str[epmelblimNR] =
132 { "no", "box size", "domain decompostion", "PME grid restriction", "maximum allowed grid scaling" };
134 struct pme_load_balancing_t {
135 gmx_bool bSepPMERanks; /**< do we have separate PME ranks? */
136 gmx_bool bActive; /**< is PME tuning active? */
137 int64_t step_rel_stop; /**< stop the tuning after this value of step_rel */
138 gmx_bool bTriggerOnDLB; /**< trigger balancing only on DD DLB */
139 gmx_bool bBalance; /**< are we in the balancing phase, i.e. trying different setups? */
140 int nstage; /**< the current maximum number of stages */
142 real cut_spacing; /**< the minimum cutoff / PME grid spacing ratio */
143 real rcut_vdw; /**< Vdw cutoff (does not change) */
144 real rcut_coulomb_start; /**< Initial electrostatics cutoff */
145 real rbufOuter_coulomb; /**< the outer pairlist buffer size */
146 real rbufOuter_vdw; /**< the outer pairlist buffer size */
147 real rbufInner_coulomb; /**< the inner pairlist buffer size */
148 real rbufInner_vdw; /**< the inner pairlist buffer size */
149 matrix box_start; /**< the initial simulation box */
150 int n; /**< the count of setup as well as the allocation size */
151 pme_setup_t *setup; /**< the PME+cutoff setups */
152 int cur; /**< the inex (in setup) of the current setup */
153 int fastest; /**< index of the fastest setup up till now */
154 int lower_limit; /**< don't go below this setup index */
155 int start; /**< start of setup index range to consider in stage>0 */
156 int end; /**< end of setup index range to consider in stage>0 */
157 int elimited; /**< was the balancing limited, uses enum above */
158 int cutoff_scheme; /**< Verlet or group cut-offs */
160 int stage; /**< the current stage */
162 int cycles_n; /**< step cycle counter cummulative count */
163 double cycles_c; /**< step cycle counter cummulative cycles */
166 /* TODO The code in this file should call this getter, rather than
167 * read bActive anywhere */
168 bool pme_loadbal_is_active(const pme_load_balancing_t *pme_lb)
170 return pme_lb != nullptr && pme_lb->bActive;
173 void pme_loadbal_init(pme_load_balancing_t **pme_lb_p,
175 const gmx::MDLogger &mdlog,
176 const t_inputrec *ir,
178 const interaction_const_t *ic,
179 const NbnxnListParameters *listParams,
184 GMX_RELEASE_ASSERT(ir->cutoff_scheme != ecutsGROUP, "PME tuning is not supported with cutoff-scheme=group (because it contains bugs)");
186 pme_load_balancing_t *pme_lb;
190 // Note that we don't (yet) support PME load balancing with LJ-PME only.
191 GMX_RELEASE_ASSERT(EEL_PME(ir->coulombtype), "pme_loadbal_init called without PME electrostatics");
192 // To avoid complexity, we require a single cut-off with PME for q+LJ.
193 // This is checked by grompp, but it doesn't hurt to check again.
194 GMX_RELEASE_ASSERT(!(EEL_PME(ir->coulombtype) && EVDW_PME(ir->vdwtype) && ir->rcoulomb != ir->rvdw), "With Coulomb and LJ PME, rcoulomb should be equal to rvdw");
198 pme_lb->bSepPMERanks = !thisRankHasDuty(cr, DUTY_PME);
200 /* Initially we turn on balancing directly on based on PP/PME imbalance */
201 pme_lb->bTriggerOnDLB = FALSE;
203 /* Any number of stages >= 2 is supported */
206 pme_lb->cutoff_scheme = ir->cutoff_scheme;
208 pme_lb->rbufOuter_coulomb = listParams->rlistOuter - ic->rcoulomb;
209 pme_lb->rbufOuter_vdw = listParams->rlistOuter - ic->rvdw;
210 pme_lb->rbufInner_coulomb = listParams->rlistInner - ic->rcoulomb;
211 pme_lb->rbufInner_vdw = listParams->rlistInner - ic->rvdw;
213 /* Scale box with Ewald wall factor; note that we pmedata->boxScaler
214 * can't always usedd as it's not available with separate PME ranks.
216 EwaldBoxZScaler boxScaler(*ir);
217 boxScaler.scaleBox(box, pme_lb->box_start);
220 snew(pme_lb->setup, pme_lb->n);
222 pme_lb->rcut_vdw = ic->rvdw;
223 pme_lb->rcut_coulomb_start = ir->rcoulomb;
226 pme_lb->setup[0].rcut_coulomb = ic->rcoulomb;
227 pme_lb->setup[0].rlistOuter = listParams->rlistOuter;
228 pme_lb->setup[0].rlistInner = listParams->rlistInner;
229 pme_lb->setup[0].grid[XX] = ir->nkx;
230 pme_lb->setup[0].grid[YY] = ir->nky;
231 pme_lb->setup[0].grid[ZZ] = ir->nkz;
232 pme_lb->setup[0].ewaldcoeff_q = ic->ewaldcoeff_q;
233 pme_lb->setup[0].ewaldcoeff_lj = ic->ewaldcoeff_lj;
235 if (!pme_lb->bSepPMERanks)
237 GMX_RELEASE_ASSERT(pmedata, "On ranks doing both PP and PME we need a valid pmedata object");
238 pme_lb->setup[0].pmedata = pmedata;
242 for (d = 0; d < DIM; d++)
244 sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
250 pme_lb->setup[0].spacing = spm;
252 if (ir->fourier_spacing > 0)
254 pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
258 pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
264 pme_lb->lower_limit = 0;
267 pme_lb->elimited = epmelblimNO;
269 pme_lb->cycles_n = 0;
270 pme_lb->cycles_c = 0;
272 if (!wallcycle_have_counter())
274 GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use PME-PP balancing.");
277 /* Tune with GPUs and/or separate PME ranks.
278 * When running only on a CPU without PME ranks, PME tuning will only help
279 * with small numbers of atoms in the cut-off sphere.
281 pme_lb->bActive = (wallcycle_have_counter() && (bUseGPU ||
282 pme_lb->bSepPMERanks));
284 /* With GPUs and no separate PME ranks we can't measure the PP/PME
285 * imbalance, so we start balancing right away.
286 * Otherwise we only start balancing after we observe imbalance.
288 pme_lb->bBalance = (pme_lb->bActive && (bUseGPU && !pme_lb->bSepPMERanks));
290 pme_lb->step_rel_stop = PMETunePeriod*ir->nstlist;
292 /* Delay DD load balancing when GPUs are used */
293 if (pme_lb->bActive && DOMAINDECOMP(cr) && cr->dd->nnodes > 1 && bUseGPU)
295 /* Lock DLB=auto to off (does nothing when DLB=yes/no.
296 * With GPUs and separate PME nodes, we want to first
297 * do PME tuning without DLB, since DLB might limit
298 * the cut-off, which never improves performance.
299 * We allow for DLB + PME tuning after a first round of tuning.
302 if (dd_dlb_is_locked(cr->dd))
304 GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: DLB will not turn on during the first phase of PME tuning");
310 *bPrinting = pme_lb->bBalance;
313 /*! \brief Try to increase the cutoff during load balancing */
314 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t *pme_lb,
316 const gmx_domdec_t *dd)
320 real tmpr_coulomb, tmpr_vdw;
324 /* Try to add a new setup with next larger cut-off to the list */
326 srenew(pme_lb->setup, pme_lb->n);
327 set = &pme_lb->setup[pme_lb->n-1];
328 set->pmedata = nullptr;
330 NumPmeDomains numPmeDomains = getNumPmeDomains(dd);
335 /* Avoid infinite while loop, which can occur at the minimum grid size.
336 * Note that in practice load balancing will stop before this point.
337 * The factor 2.1 allows for the extreme case in which only grids
338 * of powers of 2 are allowed (the current code supports more grids).
348 clear_ivec(set->grid);
349 sp = calcFftGrid(nullptr, pme_lb->box_start,
350 fac*pme_lb->setup[pme_lb->cur].spacing,
351 minimalPmeGridSize(pme_order),
356 /* As here we can't easily check if one of the PME ranks
357 * uses threading, we do a conservative grid check.
358 * This means we can't use pme_order or less grid lines
359 * per PME rank along x, which is not a strong restriction.
361 grid_ok = gmx_pme_check_restrictions(pme_order,
362 set->grid[XX], set->grid[YY], set->grid[ZZ],
367 while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
369 set->rcut_coulomb = pme_lb->cut_spacing*sp;
370 if (set->rcut_coulomb < pme_lb->rcut_coulomb_start)
372 /* This is unlikely, but can happen when e.g. continuing from
373 * a checkpoint after equilibration where the box shrank a lot.
374 * We want to avoid rcoulomb getting smaller than rvdw
375 * and there might be more issues with decreasing rcoulomb.
377 set->rcut_coulomb = pme_lb->rcut_coulomb_start;
380 if (pme_lb->cutoff_scheme == ecutsVERLET)
382 /* Never decrease the Coulomb and VdW list buffers */
383 set->rlistOuter = std::max(set->rcut_coulomb + pme_lb->rbufOuter_coulomb,
384 pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw);
385 set->rlistInner = std::max(set->rcut_coulomb + pme_lb->rbufInner_coulomb,
386 pme_lb->rcut_vdw + pme_lb->rbufInner_vdw);
390 /* TODO Remove these lines and pme_lb->cutoff_scheme */
391 tmpr_coulomb = set->rcut_coulomb + pme_lb->rbufOuter_coulomb;
392 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw;
393 /* Two (known) bugs with cutoff-scheme=group here:
394 * - This modification of rlist results in incorrect DD comunication.
395 * - We should set fr->bTwinRange = (fr->rlistlong > fr->rlist).
397 set->rlistOuter = std::min(tmpr_coulomb, tmpr_vdw);
398 set->rlistInner = set->rlistOuter;
402 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
403 set->grid_efficiency = 1;
404 for (d = 0; d < DIM; d++)
406 set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
408 /* The Ewald coefficient is inversly proportional to the cut-off */
410 pme_lb->setup[0].ewaldcoeff_q*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
411 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
413 pme_lb->setup[0].ewaldcoeff_lj*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
420 fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
421 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
426 /*! \brief Print the PME grid */
427 static void print_grid(FILE *fp_err, FILE *fp_log,
430 const pme_setup_t *set,
433 char buf[STRLEN], buft[STRLEN];
437 sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
443 sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
445 desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
447 if (fp_err != nullptr)
449 fprintf(fp_err, "\r%s\n", buf);
452 if (fp_log != nullptr)
454 fprintf(fp_log, "%s\n", buf);
458 /*! \brief Return the index of the last setup used in PME load balancing */
459 static int pme_loadbal_end(pme_load_balancing_t *pme_lb)
461 /* In the initial stage only n is set; end is not set yet */
472 /*! \brief Print descriptive string about what limits PME load balancing */
473 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
475 pme_load_balancing_t *pme_lb)
477 char buf[STRLEN], sbuf[22];
479 sprintf(buf, "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
480 gmx_step_str(step, sbuf),
481 pmelblim_str[pme_lb->elimited],
482 pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
483 if (fp_err != nullptr)
485 fprintf(fp_err, "\r%s\n", buf);
488 if (fp_log != nullptr)
490 fprintf(fp_log, "%s\n", buf);
494 /*! \brief Switch load balancing to stage 1
496 * In this stage, only reasonably fast setups are run again. */
497 static void switch_to_stage1(pme_load_balancing_t *pme_lb)
499 /* Increase start until we find a setup that is not slower than
500 * maxRelativeSlowdownAccepted times the fastest setup.
502 pme_lb->start = pme_lb->lower_limit;
503 while (pme_lb->start + 1 < pme_lb->n &&
504 (pme_lb->setup[pme_lb->start].count == 0 ||
505 pme_lb->setup[pme_lb->start].cycles >
506 pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted))
510 /* While increasing start, we might have skipped setups that we did not
511 * time during stage 0. We want to extend the range for stage 1 to include
512 * any skipped setups that lie between setups that were measured to be
513 * acceptably fast and too slow.
515 while (pme_lb->start > pme_lb->lower_limit &&
516 pme_lb->setup[pme_lb->start - 1].count == 0)
521 /* Decrease end only with setups that we timed and that are slow. */
522 pme_lb->end = pme_lb->n;
523 if (pme_lb->setup[pme_lb->end - 1].count > 0 &&
524 pme_lb->setup[pme_lb->end - 1].cycles >
525 pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted)
532 /* Next we want to choose setup pme_lb->end-1, but as we will decrease
533 * pme_lb->cur by one right after returning, we set cur to end.
535 pme_lb->cur = pme_lb->end;
538 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
540 * The adjustment is done to generate a different non-bonded PP and PME load.
541 * With separate PME ranks (PP and PME on different processes) or with
542 * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
543 * and changing the load will affect the load balance and performance.
544 * The total time for a set of integration steps is monitored and a range
545 * of grid/cut-off setups is scanned. After calling pme_load_balance many
546 * times and acquiring enough statistics, the best performing setup is chosen.
547 * Here we try to take into account fluctuations and changes due to external
548 * factors as well as DD load balancing.
551 pme_load_balance(pme_load_balancing_t *pme_lb,
555 const gmx::MDLogger &mdlog,
556 const t_inputrec *ir,
559 interaction_const_t *ic,
560 struct nonbonded_verlet_t *nbv,
561 struct gmx_pme_t ** pmedata,
567 char buf[STRLEN], sbuf[22];
572 gmx_sumd(1, &cycles, cr);
573 cycles /= cr->nnodes;
576 set = &pme_lb->setup[pme_lb->cur];
579 rtab = ir->rlist + ir->tabext;
581 if (set->count % 2 == 1)
583 /* Skip the first cycle, because the first step after a switch
584 * is much slower due to allocation and/or caching effects.
589 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
590 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
594 set->cycles = cycles;
598 if (cycles*maxFluctuationAccepted < set->cycles &&
599 pme_lb->stage == pme_lb->nstage - 1)
601 /* The performance went up a lot (due to e.g. DD load balancing).
602 * Add a stage, keep the minima, but rescan all setups.
608 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
609 "Increased the number stages to %d"
610 " and ignoring the previous performance\n",
611 set->grid[XX], set->grid[YY], set->grid[ZZ],
612 set->cycles*1e-6, cycles*1e-6, maxFluctuationAccepted,
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->n = 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 < pme_lb->n)
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->ePBC != epbcNONE)
686 OK = (gmx::square(pme_lb->setup[pme_lb->cur+1].rlistOuter)
687 <= max_cutoff2(ir->ePBC, state->box));
690 pme_lb->elimited = epmelblimBOX;
698 if (DOMAINDECOMP(cr))
700 OK = change_dd_cutoff(cr, state, ir,
701 pme_lb->setup[pme_lb->cur].rlistOuter);
704 /* Failed: do not use this setup */
706 pme_lb->elimited = epmelblimDD;
712 /* We hit the upper limit for the cut-off,
713 * the setup should not go further than cur.
715 pme_lb->n = pme_lb->cur + 1;
716 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
717 /* Switch to the next stage */
718 switch_to_stage1(pme_lb);
722 !(pme_lb->setup[pme_lb->cur].grid[XX]*
723 pme_lb->setup[pme_lb->cur].grid[YY]*
724 pme_lb->setup[pme_lb->cur].grid[ZZ] <
725 gridsize_start*gridpointsScaleFactor
727 pme_lb->setup[pme_lb->cur].grid_efficiency <
728 pme_lb->setup[pme_lb->cur-1].grid_efficiency*relativeEfficiencyFactor));
731 if (pme_lb->stage > 0 && pme_lb->end == 1)
733 pme_lb->cur = pme_lb->lower_limit;
734 pme_lb->stage = pme_lb->nstage;
736 else if (pme_lb->stage > 0 && pme_lb->end > 1)
738 /* If stage = nstage-1:
739 * scan over all setups, rerunning only those setups
740 * which are not much slower than the fastest
743 * Note that we loop backward to minimize the risk of the cut-off
744 * getting limited by DD DLB, since the DLB cut-off limit is set
745 * to the fastest PME setup.
749 if (pme_lb->cur > pme_lb->start)
757 pme_lb->cur = pme_lb->end - 1;
760 while (pme_lb->stage == pme_lb->nstage - 1 &&
761 pme_lb->setup[pme_lb->cur].count > 0 &&
762 pme_lb->setup[pme_lb->cur].cycles > cycles_fast*maxRelativeSlowdownAccepted);
764 if (pme_lb->stage == pme_lb->nstage)
766 /* We are done optimizing, use the fastest setup we found */
767 pme_lb->cur = pme_lb->fastest;
771 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
773 OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlistOuter);
776 /* For some reason the chosen cut-off is incompatible with DD.
777 * We should continue scanning a more limited range of cut-off's.
779 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
781 /* stage=nstage says we're finished, but we should continue
782 * balancing, so we set back stage which was just incremented.
786 if (pme_lb->cur <= pme_lb->fastest)
788 /* This should not happen, as we set limits on the DLB bounds.
789 * But we implement a complete failsafe solution anyhow.
791 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
792 "The fastest PP/PME load balancing setting (cutoff %.3d nm) is no longer available due to DD DLB or box size limitations", pme_lb->fastest);
793 pme_lb->fastest = pme_lb->lower_limit;
794 pme_lb->start = pme_lb->lower_limit;
796 /* Limit the range to below the current cut-off, scan from start */
797 pme_lb->end = pme_lb->cur;
798 pme_lb->cur = pme_lb->start;
799 pme_lb->elimited = epmelblimDD;
800 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
804 /* Change the Coulomb cut-off and the PME grid */
806 set = &pme_lb->setup[pme_lb->cur];
808 NbnxnListParameters *listParams = nbv->listParams.get();
810 ic->rcoulomb = set->rcut_coulomb;
811 listParams->rlistOuter = set->rlistOuter;
812 listParams->rlistInner = set->rlistInner;
813 ic->ewaldcoeff_q = set->ewaldcoeff_q;
814 /* TODO: centralize the code that sets the potentials shifts */
815 if (ic->coulomb_modifier == eintmodPOTSHIFT)
817 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
818 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
820 if (EVDW_PME(ic->vdwtype))
822 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
823 ic->rvdw = set->rcut_coulomb;
824 ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
825 if (ic->vdw_modifier == eintmodPOTSHIFT)
829 ic->dispersion_shift.cpot = -1.0/gmx::power6(static_cast<double>(ic->rvdw));
830 ic->repulsion_shift.cpot = -1.0/gmx::power12(static_cast<double>(ic->rvdw));
831 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
832 crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
833 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
837 /* We always re-initialize the tables whether they are used or not */
838 init_interaction_const_tables(nullptr, ic, rtab);
840 nbnxn_gpu_pme_loadbal_update_param(nbv, ic, listParams);
842 if (!pme_lb->bSepPMERanks)
845 * CPU PME keeps a list of allocated pmedata's, that's why pme_lb->setup[pme_lb->cur].pmedata is not always nullptr.
846 * GPU PME, however, currently needs the gmx_pme_reinit always called on load balancing
847 * (pme_gpu_reinit might be not sufficiently decoupled from gmx_pme_init).
848 * This can lead to a lot of reallocations for PME GPU.
849 * Would be nicer if the allocated grid list was hidden within a single pmedata structure.
851 if ((pme_lb->setup[pme_lb->cur].pmedata == nullptr) || pme_gpu_task_enabled(pme_lb->setup[pme_lb->cur].pmedata))
853 /* Generate a new PME data structure,
854 * copying part of the old pointers.
856 gmx_pme_reinit(&set->pmedata,
857 cr, pme_lb->setup[0].pmedata, ir,
858 set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
860 *pmedata = set->pmedata;
864 /* Tell our PME-only rank to switch grid */
865 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
870 print_grid(nullptr, debug, "", "switched to", set, -1);
873 if (pme_lb->stage == pme_lb->nstage)
875 print_grid(fp_err, fp_log, "", "optimal", set, -1);
879 /*! \brief Prepare for another round of PME load balancing
881 * \param[in,out] pme_lb Pointer to PME load balancing struct
882 * \param[in] bDlbUnlocked TRUE is DLB was locked and is now unlocked
884 * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
885 * the PP/PME balance might change and re-balancing can improve performance.
886 * This function adds 2 stages and adjusts the considered setup range.
888 static void continue_pme_loadbal(pme_load_balancing_t *pme_lb,
889 gmx_bool bDlbUnlocked)
891 /* Add 2 tuning stages, keep the detected end of the setup range */
893 if (bDlbUnlocked && pme_lb->bSepPMERanks)
895 /* With separate PME ranks, DLB should always lower the PP load and
896 * can only increase the PME load (more communication and imbalance),
897 * so we only need to scan longer cut-off's.
899 pme_lb->lower_limit = pme_lb->cur;
901 pme_lb->start = pme_lb->lower_limit;
904 void pme_loadbal_do(pme_load_balancing_t *pme_lb,
908 const gmx::MDLogger &mdlog,
909 const t_inputrec *ir,
912 gmx_wallcycle_t wcycle,
920 assert(pme_lb != nullptr);
922 if (!pme_lb->bActive)
927 n_prev = pme_lb->cycles_n;
928 cycles_prev = pme_lb->cycles_c;
929 wallcycle_get(wcycle, ewcSTEP, &pme_lb->cycles_n, &pme_lb->cycles_c);
931 /* Before the first step we haven't done any steps yet.
932 * Also handle cases where ir->init_step % ir->nstlist != 0.
934 if (pme_lb->cycles_n < ir->nstlist)
938 /* Sanity check, we expect nstlist cycle counts */
939 if (pme_lb->cycles_n - n_prev != ir->nstlist)
941 /* We could return here, but it's safer to issue an error and quit */
942 gmx_incons("pme_loadbal_do called at an interval != nstlist");
945 /* PME grid + cut-off optimization with GPUs or PME ranks */
946 if (!pme_lb->bBalance && pme_lb->bSepPMERanks)
948 if (pme_lb->bTriggerOnDLB)
950 pme_lb->bBalance = dd_dlb_is_on(cr->dd);
952 /* We should ignore the first timing to avoid timing allocation
953 * overhead. And since the PME load balancing is called just
954 * before DD repartitioning, the ratio returned by dd_pme_f_ratio
955 * is not over the last nstlist steps, but the nstlist steps before
956 * that. So the first useful ratio is available at step_rel=3*nstlist.
958 else if (step_rel >= 3*ir->nstlist)
960 if (DDMASTER(cr->dd))
962 /* If PME rank load is too high, start tuning */
964 (dd_pme_f_ratio(cr->dd) >= loadBalanceTriggerFactor);
966 dd_bcast(cr->dd, sizeof(gmx_bool), &pme_lb->bBalance);
969 pme_lb->bActive = (pme_lb->bBalance ||
970 step_rel <= pme_lb->step_rel_stop);
973 /* The location in the code of this balancing termination is strange.
974 * You would expect to have it after the call to pme_load_balance()
975 * below, since there pme_lb->stage is updated.
976 * But when terminating directly after deciding on and selecting the
977 * optimal setup, DLB will turn on right away if it was locked before.
978 * This might be due to PME reinitialization. So we check stage here
979 * to allow for another nstlist steps with DLB locked to stabilize
982 if (pme_lb->bBalance && pme_lb->stage == pme_lb->nstage)
984 pme_lb->bBalance = FALSE;
986 if (DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
988 /* Unlock the DLB=auto, DLB is allowed to activate */
989 dd_dlb_unlock(cr->dd);
990 GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: DLB can now turn on, when beneficial");
992 /* We don't deactivate the tuning yet, since we will balance again
993 * after DLB gets turned on, if it does within PMETune_period.
995 continue_pme_loadbal(pme_lb, TRUE);
996 pme_lb->bTriggerOnDLB = TRUE;
997 pme_lb->step_rel_stop = step_rel + PMETunePeriod*ir->nstlist;
1001 /* We're completely done with PME tuning */
1002 pme_lb->bActive = FALSE;
1005 if (DOMAINDECOMP(cr))
1007 /* Set the cut-off limit to the final selected cut-off,
1008 * so we don't have artificial DLB limits.
1009 * This also ensures that we won't disable the currently
1010 * optimal setting during a second round of PME balancing.
1012 set_dd_dlb_max_cutoff(cr, fr->nbv->listParams->rlistOuter);
1016 if (pme_lb->bBalance)
1018 /* We might not have collected nstlist steps in cycles yet,
1019 * since init_step might not be a multiple of nstlist,
1020 * but the first data collected is skipped anyhow.
1022 pme_load_balance(pme_lb, cr,
1023 fp_err, fp_log, mdlog,
1024 ir, state, pme_lb->cycles_c - cycles_prev,
1025 fr->ic, fr->nbv, &fr->pmedata,
1028 /* Update deprecated rlist in forcerec to stay in sync with fr->nbv */
1029 fr->rlist = fr->nbv->listParams->rlistOuter;
1031 if (ir->eDispCorr != edispcNO)
1033 calc_enervirdiff(nullptr, ir->eDispCorr, fr);
1037 if (!pme_lb->bBalance &&
1038 (!pme_lb->bSepPMERanks || step_rel > pme_lb->step_rel_stop))
1040 /* We have just deactivated the balancing and we're not measuring PP/PME
1041 * imbalance during the first steps of the run: deactivate the tuning.
1043 pme_lb->bActive = FALSE;
1046 if (!(pme_lb->bActive) && DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
1048 /* Make sure DLB is allowed when we deactivate PME tuning */
1049 dd_dlb_unlock(cr->dd);
1050 GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: DLB can now turn on, when beneficial");
1053 *bPrinting = pme_lb->bBalance;
1056 /*! \brief Return product of the number of PME grid points in each dimension */
1057 static int pme_grid_points(const pme_setup_t *setup)
1059 return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
1062 /*! \brief Print one load-balancing setting */
1063 static void print_pme_loadbal_setting(FILE *fplog,
1065 const pme_setup_t *setup)
1068 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
1070 setup->rcut_coulomb, setup->rlistInner,
1071 setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
1072 setup->spacing, 1/setup->ewaldcoeff_q);
1075 /*! \brief Print all load-balancing settings */
1076 static void print_pme_loadbal_settings(pme_load_balancing_t *pme_lb,
1078 const gmx::MDLogger &mdlog,
1079 gmx_bool bNonBondedOnGPU)
1081 double pp_ratio, grid_ratio;
1082 real pp_ratio_temporary;
1084 pp_ratio_temporary = pme_lb->setup[pme_lb->cur].rlistInner / pme_lb->setup[0].rlistInner;
1085 pp_ratio = gmx::power3(pp_ratio_temporary);
1086 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
1087 (double)pme_grid_points(&pme_lb->setup[0]);
1089 fprintf(fplog, "\n");
1090 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
1091 fprintf(fplog, "\n");
1092 /* Here we only warn when the optimal setting is the last one */
1093 if (pme_lb->elimited != epmelblimNO &&
1094 pme_lb->cur == pme_loadbal_end(pme_lb)-1)
1096 fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
1097 pmelblim_str[pme_lb->elimited]);
1098 fprintf(fplog, " you might not have reached a good load balance.\n");
1099 if (pme_lb->elimited == epmelblimDD)
1101 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
1103 fprintf(fplog, "\n");
1105 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
1106 fprintf(fplog, " particle-particle PME\n");
1107 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
1108 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
1109 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
1110 fprintf(fplog, " cost-ratio %4.2f %4.2f\n",
1111 pp_ratio, grid_ratio);
1112 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
1114 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
1116 GMX_LOG(mdlog.warning).asParagraph().appendText(
1117 "NOTE: PME load balancing increased the non-bonded workload by more than 50%.\n"
1118 " For better performance, use (more) PME ranks (mdrun -npme),\n"
1119 " or if you are beyond the scaling limit, use fewer total ranks (or nodes).");
1123 fprintf(fplog, "\n");
1127 void pme_loadbal_done(pme_load_balancing_t *pme_lb,
1129 const gmx::MDLogger &mdlog,
1130 gmx_bool bNonBondedOnGPU)
1132 if (fplog != nullptr && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
1134 print_pme_loadbal_settings(pme_lb, fplog, mdlog, bNonBondedOnGPU);
1137 /* TODO: Here we should free all pointers in pme_lb,
1138 * but as it contains pme data structures,
1139 * we need to first make pme.c free all data.