Fix cmake policy warning
[gromacs.git] / src / gromacs / ewald / pme-load-balancing.cpp
bloba975b49788c517b1c8bba48e498109ee63f39396
1 /*
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 * 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.
35 /*! \internal \file
37 * \brief This file contains function definitions necessary for
38 * managing automatic load balance of PME calculations (Coulomb and
39 * LJ).
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_ewald
44 #include "gmxpre.h"
46 #include "pme-load-balancing.h"
48 #include <cassert>
49 #include <cmath>
51 #include <algorithm>
53 #include "gromacs/domdec/dlb.h"
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_network.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/domdec/partition.h"
58 #include "gromacs/ewald/ewald-utils.h"
59 #include "gromacs/ewald/pme.h"
60 #include "gromacs/fft/calcgrid.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/math/functions.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdlib/forcerec.h"
65 #include "gromacs/mdlib/nb_verlet.h"
66 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
67 #include "gromacs/mdlib/nbnxn_pairlist.h"
68 #include "gromacs/mdlib/sim_util.h"
69 #include "gromacs/mdtypes/commrec.h"
70 #include "gromacs/mdtypes/inputrec.h"
71 #include "gromacs/mdtypes/md_enums.h"
72 #include "gromacs/mdtypes/state.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/logger.h"
79 #include "gromacs/utility/smalloc.h"
80 #include "gromacs/utility/strconvert.h"
82 #include "pme-internal.h"
84 /*! \brief Parameters and settings for one PP-PME setup */
85 struct pme_setup_t {
86 real rcut_coulomb; /**< Coulomb cut-off */
87 real rlistOuter; /**< cut-off for the outer pair-list */
88 real rlistInner; /**< cut-off for the inner pair-list */
89 real spacing; /**< (largest) PME grid spacing */
90 ivec grid; /**< the PME grid dimensions */
91 real grid_efficiency; /**< ineffiency factor for non-uniform grids <= 1 */
92 real ewaldcoeff_q; /**< Electrostatic Ewald coefficient */
93 real ewaldcoeff_lj; /**< LJ Ewald coefficient, only for the call to send_switchgrid */
94 struct gmx_pme_t *pmedata; /**< the data structure used in the PME code */
95 int count; /**< number of times this setup has been timed */
96 double cycles; /**< the fastest time for this setup in cycles */
99 /*! \brief After 50 nstlist periods of not observing imbalance: never tune PME */
100 const int PMETunePeriod = 50;
101 /*! \brief Trigger PME load balancing at more than 5% PME overload */
102 const real loadBalanceTriggerFactor = 1.05;
103 /*! \brief Scale the grid by a most at factor 1.7.
105 * This still leaves room for about 4-4.5x decrease in grid spacing while limiting the cases where
106 * large imbalance leads to extreme cutoff scaling for marginal benefits.
108 * This should help to avoid:
109 * - large increase in power consumption for little performance gain
110 * - increasing communication volume
111 * - limiting DLB
113 const real c_maxSpacingScaling = 1.7;
114 /*! \brief In the initial scan, step by grids that are at least a factor 0.8 coarser */
115 const real gridpointsScaleFactor = 0.8;
116 /*! \brief In the initial scan, try to skip grids with uneven x/y/z spacing,
117 * checking if the "efficiency" is more than 5% worse than the previous grid.
119 const real relativeEfficiencyFactor = 1.05;
120 /*! \brief Rerun until a run is 12% slower setups than the fastest run so far */
121 const real maxRelativeSlowdownAccepted = 1.12;
122 /*! \brief If setups get more than 2% faster, do another round to avoid
123 * choosing a slower setup due to acceleration or fluctuations.
125 const real maxFluctuationAccepted = 1.02;
127 /*! \brief Enumeration whose values describe the effect limiting the load balancing */
128 enum epmelb {
129 epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimPMEGRID, epmelblimMAXSCALING, epmelblimNR
132 /*! \brief Descriptive strings matching ::epmelb */
133 static const char *pmelblim_str[epmelblimNR] =
134 { "no", "box size", "domain decompostion", "PME grid restriction", "maximum allowed grid scaling" };
136 struct pme_load_balancing_t {
137 gmx_bool bSepPMERanks; /**< do we have separate PME ranks? */
138 gmx_bool bActive; /**< is PME tuning active? */
139 int64_t step_rel_stop; /**< stop the tuning after this value of step_rel */
140 gmx_bool bTriggerOnDLB; /**< trigger balancing only on DD DLB */
141 gmx_bool bBalance; /**< are we in the balancing phase, i.e. trying different setups? */
142 int nstage; /**< the current maximum number of stages */
144 real cut_spacing; /**< the minimum cutoff / PME grid spacing ratio */
145 real rcut_vdw; /**< Vdw cutoff (does not change) */
146 real rcut_coulomb_start; /**< Initial electrostatics cutoff */
147 real rbufOuter_coulomb; /**< the outer pairlist buffer size */
148 real rbufOuter_vdw; /**< the outer pairlist buffer size */
149 real rbufInner_coulomb; /**< the inner pairlist buffer size */
150 real rbufInner_vdw; /**< the inner pairlist buffer size */
151 matrix box_start; /**< the initial simulation box */
152 std::vector<pme_setup_t> setup; /**< the PME+cutoff setups */
153 int cur; /**< the index (in setup) of the current setup */
154 int fastest; /**< index of the fastest setup up till now */
155 int lower_limit; /**< don't go below this setup index */
156 int start; /**< start of setup index range to consider in stage>0 */
157 int end; /**< end of setup index range to consider in stage>0 */
158 int elimited; /**< was the balancing limited, uses enum above */
159 int cutoff_scheme; /**< Verlet or group cut-offs */
161 int stage; /**< the current stage */
163 int cycles_n; /**< step cycle counter cummulative count */
164 double cycles_c; /**< step cycle counter cummulative cycles */
167 /* TODO The code in this file should call this getter, rather than
168 * read bActive anywhere */
169 bool pme_loadbal_is_active(const pme_load_balancing_t *pme_lb)
171 return pme_lb != nullptr && pme_lb->bActive;
174 // TODO Return a unique_ptr to pme_load_balancing_t
175 void pme_loadbal_init(pme_load_balancing_t **pme_lb_p,
176 t_commrec *cr,
177 const gmx::MDLogger &mdlog,
178 const t_inputrec &ir,
179 const matrix box,
180 const interaction_const_t &ic,
181 const NbnxnListParameters &listParams,
182 gmx_pme_t *pmedata,
183 gmx_bool bUseGPU,
184 gmx_bool *bPrinting)
186 GMX_RELEASE_ASSERT(ir.cutoff_scheme != ecutsGROUP, "PME tuning is not supported with cutoff-scheme=group (because it contains bugs)");
188 pme_load_balancing_t *pme_lb;
189 real spm, sp;
190 int d;
192 // Note that we don't (yet) support PME load balancing with LJ-PME only.
193 GMX_RELEASE_ASSERT(EEL_PME(ir.coulombtype), "pme_loadbal_init called without PME electrostatics");
194 // To avoid complexity, we require a single cut-off with PME for q+LJ.
195 // This is checked by grompp, but it doesn't hurt to check again.
196 GMX_RELEASE_ASSERT(!(EEL_PME(ir.coulombtype) && EVDW_PME(ir.vdwtype) && ir.rcoulomb != ir.rvdw), "With Coulomb and LJ PME, rcoulomb should be equal to rvdw");
198 pme_lb = new pme_load_balancing_t;
200 pme_lb->bSepPMERanks = !thisRankHasDuty(cr, DUTY_PME);
202 /* Initially we turn on balancing directly on based on PP/PME imbalance */
203 pme_lb->bTriggerOnDLB = FALSE;
205 /* Any number of stages >= 2 is supported */
206 pme_lb->nstage = 2;
208 pme_lb->cutoff_scheme = ir.cutoff_scheme;
210 pme_lb->rbufOuter_coulomb = listParams.rlistOuter - ic.rcoulomb;
211 pme_lb->rbufOuter_vdw = listParams.rlistOuter - ic.rvdw;
212 pme_lb->rbufInner_coulomb = listParams.rlistInner - ic.rcoulomb;
213 pme_lb->rbufInner_vdw = listParams.rlistInner - ic.rvdw;
215 /* Scale box with Ewald wall factor; note that we pmedata->boxScaler
216 * can't always usedd as it's not available with separate PME ranks.
218 EwaldBoxZScaler boxScaler(ir);
219 boxScaler.scaleBox(box, pme_lb->box_start);
221 pme_lb->setup.resize(1);
223 pme_lb->rcut_vdw = ic.rvdw;
224 pme_lb->rcut_coulomb_start = ir.rcoulomb;
226 pme_lb->cur = 0;
227 pme_lb->setup[0].rcut_coulomb = ic.rcoulomb;
228 pme_lb->setup[0].rlistOuter = listParams.rlistOuter;
229 pme_lb->setup[0].rlistInner = listParams.rlistInner;
230 pme_lb->setup[0].grid[XX] = ir.nkx;
231 pme_lb->setup[0].grid[YY] = ir.nky;
232 pme_lb->setup[0].grid[ZZ] = ir.nkz;
233 pme_lb->setup[0].ewaldcoeff_q = ic.ewaldcoeff_q;
234 pme_lb->setup[0].ewaldcoeff_lj = ic.ewaldcoeff_lj;
236 if (!pme_lb->bSepPMERanks)
238 GMX_RELEASE_ASSERT(pmedata, "On ranks doing both PP and PME we need a valid pmedata object");
239 pme_lb->setup[0].pmedata = pmedata;
242 spm = 0;
243 for (d = 0; d < DIM; d++)
245 sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
246 if (sp > spm)
248 spm = sp;
251 pme_lb->setup[0].spacing = spm;
253 if (ir.fourier_spacing > 0)
255 pme_lb->cut_spacing = ir.rcoulomb/ir.fourier_spacing;
257 else
259 pme_lb->cut_spacing = ir.rcoulomb/pme_lb->setup[0].spacing;
262 pme_lb->stage = 0;
264 pme_lb->fastest = 0;
265 pme_lb->lower_limit = 0;
266 pme_lb->start = 0;
267 pme_lb->end = 0;
268 pme_lb->elimited = epmelblimNO;
270 pme_lb->cycles_n = 0;
271 pme_lb->cycles_c = 0;
273 if (!wallcycle_have_counter())
275 GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use PME-PP balancing.");
278 /* Tune with GPUs and/or separate PME ranks.
279 * When running only on a CPU without PME ranks, PME tuning will only help
280 * with small numbers of atoms in the cut-off sphere.
282 pme_lb->bActive = (wallcycle_have_counter() && (bUseGPU ||
283 pme_lb->bSepPMERanks));
285 /* With GPUs and no separate PME ranks we can't measure the PP/PME
286 * imbalance, so we start balancing right away.
287 * Otherwise we only start balancing after we observe imbalance.
289 pme_lb->bBalance = (pme_lb->bActive && (bUseGPU && !pme_lb->bSepPMERanks));
291 pme_lb->step_rel_stop = PMETunePeriod*ir.nstlist;
293 /* Delay DD load balancing when GPUs are used */
294 if (pme_lb->bActive && DOMAINDECOMP(cr) && cr->dd->nnodes > 1 && bUseGPU)
296 /* Lock DLB=auto to off (does nothing when DLB=yes/no.
297 * With GPUs and separate PME nodes, we want to first
298 * do PME tuning without DLB, since DLB might limit
299 * the cut-off, which never improves performance.
300 * We allow for DLB + PME tuning after a first round of tuning.
302 dd_dlb_lock(cr->dd);
303 if (dd_dlb_is_locked(cr->dd))
305 GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: DLB will not turn on during the first phase of PME tuning");
309 *pme_lb_p = pme_lb;
311 *bPrinting = pme_lb->bBalance;
314 /*! \brief Try to increase the cutoff during load balancing */
315 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t *pme_lb,
316 int pme_order,
317 const gmx_domdec_t *dd)
319 real fac, sp;
320 real tmpr_coulomb, tmpr_vdw;
321 int d;
322 bool grid_ok;
324 /* Try to add a new setup with next larger cut-off to the list */
325 pme_setup_t set;
327 set.pmedata = nullptr;
329 NumPmeDomains numPmeDomains = getNumPmeDomains(dd);
331 fac = 1;
334 /* Avoid infinite while loop, which can occur at the minimum grid size.
335 * Note that in practice load balancing will stop before this point.
336 * The factor 2.1 allows for the extreme case in which only grids
337 * of powers of 2 are allowed (the current code supports more grids).
339 if (fac > 2.1)
341 return FALSE;
344 fac *= 1.01;
345 clear_ivec(set.grid);
346 sp = calcFftGrid(nullptr, pme_lb->box_start,
347 fac*pme_lb->setup[pme_lb->cur].spacing,
348 minimalPmeGridSize(pme_order),
349 &set.grid[XX],
350 &set.grid[YY],
351 &set.grid[ZZ]);
353 /* As here we can't easily check if one of the PME ranks
354 * uses threading, we do a conservative grid check.
355 * This means we can't use pme_order or less grid lines
356 * per PME rank along x, which is not a strong restriction.
358 grid_ok = gmx_pme_check_restrictions(pme_order,
359 set.grid[XX], set.grid[YY], set.grid[ZZ],
360 numPmeDomains.x,
361 true,
362 false);
364 while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
366 set.rcut_coulomb = pme_lb->cut_spacing*sp;
367 if (set.rcut_coulomb < pme_lb->rcut_coulomb_start)
369 /* This is unlikely, but can happen when e.g. continuing from
370 * a checkpoint after equilibration where the box shrank a lot.
371 * We want to avoid rcoulomb getting smaller than rvdw
372 * and there might be more issues with decreasing rcoulomb.
374 set.rcut_coulomb = pme_lb->rcut_coulomb_start;
377 if (pme_lb->cutoff_scheme == ecutsVERLET)
379 /* Never decrease the Coulomb and VdW list buffers */
380 set.rlistOuter = std::max(set.rcut_coulomb + pme_lb->rbufOuter_coulomb,
381 pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw);
382 set.rlistInner = std::max(set.rcut_coulomb + pme_lb->rbufInner_coulomb,
383 pme_lb->rcut_vdw + pme_lb->rbufInner_vdw);
385 else
387 /* TODO Remove these lines and pme_lb->cutoff_scheme */
388 tmpr_coulomb = set.rcut_coulomb + pme_lb->rbufOuter_coulomb;
389 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw;
390 /* Two (known) bugs with cutoff-scheme=group here:
391 * - This modification of rlist results in incorrect DD comunication.
392 * - We should set fr->bTwinRange = (fr->rlistlong > fr->rlist).
394 set.rlistOuter = std::min(tmpr_coulomb, tmpr_vdw);
395 set.rlistInner = set.rlistOuter;
398 set.spacing = sp;
399 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
400 set.grid_efficiency = 1;
401 for (d = 0; d < DIM; d++)
403 set.grid_efficiency *= (set.grid[d]*sp)/norm(pme_lb->box_start[d]);
405 /* The Ewald coefficient is inversly proportional to the cut-off */
406 set.ewaldcoeff_q =
407 pme_lb->setup[0].ewaldcoeff_q*pme_lb->setup[0].rcut_coulomb/set.rcut_coulomb;
408 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
409 set.ewaldcoeff_lj =
410 pme_lb->setup[0].ewaldcoeff_lj*pme_lb->setup[0].rcut_coulomb/set.rcut_coulomb;
412 set.count = 0;
413 set.cycles = 0;
415 if (debug)
417 fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
418 set.grid[XX], set.grid[YY], set.grid[ZZ], set.rcut_coulomb);
420 pme_lb->setup.push_back(set);
421 return TRUE;
424 /*! \brief Print the PME grid */
425 static void print_grid(FILE *fp_err, FILE *fp_log,
426 const char *pre,
427 const char *desc,
428 const pme_setup_t *set,
429 double cycles)
431 auto buf = gmx::formatString("%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f",
432 pre, desc,
433 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
434 if (cycles >= 0)
436 buf += gmx::formatString(": %.1f M-cycles", cycles*1e-6);
438 if (fp_err != nullptr)
440 fprintf(fp_err, "\r%s\n", buf.c_str());
441 fflush(fp_err);
443 if (fp_log != nullptr)
445 fprintf(fp_log, "%s\n", buf.c_str());
449 /*! \brief Return the index of the last setup used in PME load balancing */
450 static int pme_loadbal_end(pme_load_balancing_t *pme_lb)
452 /* In the initial stage only n is set; end is not set yet */
453 if (pme_lb->end > 0)
455 return pme_lb->end;
457 else
459 return pme_lb->setup.size();
463 /*! \brief Print descriptive string about what limits PME load balancing */
464 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
465 int64_t step,
466 pme_load_balancing_t *pme_lb)
468 auto buf = gmx::formatString("step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
469 gmx::int64ToString(step).c_str(),
470 pmelblim_str[pme_lb->elimited],
471 pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
472 if (fp_err != nullptr)
474 fprintf(fp_err, "\r%s\n", buf.c_str());
475 fflush(fp_err);
477 if (fp_log != nullptr)
479 fprintf(fp_log, "%s\n", buf.c_str());
483 /*! \brief Switch load balancing to stage 1
485 * In this stage, only reasonably fast setups are run again. */
486 static void switch_to_stage1(pme_load_balancing_t *pme_lb)
488 /* Increase start until we find a setup that is not slower than
489 * maxRelativeSlowdownAccepted times the fastest setup.
491 pme_lb->start = pme_lb->lower_limit;
492 while (pme_lb->start + 1 < static_cast<int>(pme_lb->setup.size()) &&
493 (pme_lb->setup[pme_lb->start].count == 0 ||
494 pme_lb->setup[pme_lb->start].cycles >
495 pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted))
497 pme_lb->start++;
499 /* While increasing start, we might have skipped setups that we did not
500 * time during stage 0. We want to extend the range for stage 1 to include
501 * any skipped setups that lie between setups that were measured to be
502 * acceptably fast and too slow.
504 while (pme_lb->start > pme_lb->lower_limit &&
505 pme_lb->setup[pme_lb->start - 1].count == 0)
507 pme_lb->start--;
510 /* Decrease end only with setups that we timed and that are slow. */
511 pme_lb->end = pme_lb->setup.size();
512 if (pme_lb->setup[pme_lb->end - 1].count > 0 &&
513 pme_lb->setup[pme_lb->end - 1].cycles >
514 pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted)
516 pme_lb->end--;
519 pme_lb->stage = 1;
521 /* Next we want to choose setup pme_lb->end-1, but as we will decrease
522 * pme_lb->cur by one right after returning, we set cur to end.
524 pme_lb->cur = pme_lb->end;
527 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
529 * The adjustment is done to generate a different non-bonded PP and PME load.
530 * With separate PME ranks (PP and PME on different processes) or with
531 * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
532 * and changing the load will affect the load balance and performance.
533 * The total time for a set of integration steps is monitored and a range
534 * of grid/cut-off setups is scanned. After calling pme_load_balance many
535 * times and acquiring enough statistics, the best performing setup is chosen.
536 * Here we try to take into account fluctuations and changes due to external
537 * factors as well as DD load balancing.
539 static void
540 pme_load_balance(pme_load_balancing_t *pme_lb,
541 t_commrec *cr,
542 FILE *fp_err,
543 FILE *fp_log,
544 const gmx::MDLogger &mdlog,
545 const t_inputrec &ir,
546 const t_state &state,
547 double cycles,
548 interaction_const_t *ic,
549 struct nonbonded_verlet_t *nbv,
550 struct gmx_pme_t ** pmedata,
551 int64_t step)
553 gmx_bool OK;
554 pme_setup_t *set;
555 double cycles_fast;
556 char buf[STRLEN], sbuf[22];
557 real rtab;
559 if (PAR(cr))
561 gmx_sumd(1, &cycles, cr);
562 cycles /= cr->nnodes;
565 set = &pme_lb->setup[pme_lb->cur];
566 set->count++;
568 rtab = ir.rlist + ir.tabext;
570 if (set->count % 2 == 1)
572 /* Skip the first cycle, because the first step after a switch
573 * is much slower due to allocation and/or caching effects.
575 return;
578 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
579 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
581 if (set->count <= 2)
583 set->cycles = cycles;
585 else
587 if (cycles*maxFluctuationAccepted < set->cycles &&
588 pme_lb->stage == pme_lb->nstage - 1)
590 /* The performance went up a lot (due to e.g. DD load balancing).
591 * Add a stage, keep the minima, but rescan all setups.
593 pme_lb->nstage++;
595 if (debug)
597 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
598 "Increased the number stages to %d"
599 " and ignoring the previous performance\n",
600 set->grid[XX], set->grid[YY], set->grid[ZZ],
601 set->cycles*1e-6, cycles*1e-6, maxFluctuationAccepted,
602 pme_lb->nstage);
605 set->cycles = std::min(set->cycles, cycles);
608 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
610 pme_lb->fastest = pme_lb->cur;
612 if (DOMAINDECOMP(cr))
614 /* We found a new fastest setting, ensure that with subsequent
615 * shorter cut-off's the dynamic load balancing does not make
616 * the use of the current cut-off impossible. This solution is
617 * a trade-off, as the PME load balancing and DD domain size
618 * load balancing can interact in complex ways.
619 * With the Verlet kernels, DD load imbalance will usually be
620 * mainly due to bonded interaction imbalance, which will often
621 * quickly push the domain boundaries beyond the limit for the
622 * optimal, PME load balanced, cut-off. But it could be that
623 * better overal performance can be obtained with a slightly
624 * shorter cut-off and better DD load balancing.
626 set_dd_dlb_max_cutoff(cr, pme_lb->setup[pme_lb->fastest].rlistOuter);
629 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
631 /* Check in stage 0 if we should stop scanning grids.
632 * Stop when the time is more than maxRelativeSlowDownAccepted longer than the fastest.
634 if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
635 cycles > pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted)
637 pme_lb->setup.resize(pme_lb->cur + 1);
638 /* Done with scanning, go to stage 1 */
639 switch_to_stage1(pme_lb);
642 if (pme_lb->stage == 0)
644 int gridsize_start;
646 gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
650 if (pme_lb->cur+1 < static_cast<int>(pme_lb->setup.size()))
652 /* We had already generated the next setup */
653 OK = TRUE;
655 else
657 /* Find the next setup */
658 OK = pme_loadbal_increase_cutoff(pme_lb, ir.pme_order, cr->dd);
660 if (!OK)
662 pme_lb->elimited = epmelblimPMEGRID;
666 if (OK &&
667 pme_lb->setup[pme_lb->cur+1].spacing > c_maxSpacingScaling*pme_lb->setup[0].spacing)
669 OK = FALSE;
670 pme_lb->elimited = epmelblimMAXSCALING;
673 if (OK && ir.ePBC != epbcNONE)
675 OK = (gmx::square(pme_lb->setup[pme_lb->cur+1].rlistOuter)
676 <= max_cutoff2(ir.ePBC, state.box));
677 if (!OK)
679 pme_lb->elimited = epmelblimBOX;
683 if (OK)
685 pme_lb->cur++;
687 if (DOMAINDECOMP(cr))
689 OK = change_dd_cutoff(cr, state,
690 pme_lb->setup[pme_lb->cur].rlistOuter);
691 if (!OK)
693 /* Failed: do not use this setup */
694 pme_lb->cur--;
695 pme_lb->elimited = epmelblimDD;
699 if (!OK)
701 /* We hit the upper limit for the cut-off,
702 * the setup should not go further than cur.
704 pme_lb->setup.resize(pme_lb->cur + 1);
705 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
706 /* Switch to the next stage */
707 switch_to_stage1(pme_lb);
710 while (OK &&
711 !(pme_lb->setup[pme_lb->cur].grid[XX]*
712 pme_lb->setup[pme_lb->cur].grid[YY]*
713 pme_lb->setup[pme_lb->cur].grid[ZZ] <
714 gridsize_start*gridpointsScaleFactor
716 pme_lb->setup[pme_lb->cur].grid_efficiency <
717 pme_lb->setup[pme_lb->cur-1].grid_efficiency*relativeEfficiencyFactor));
720 if (pme_lb->stage > 0 && pme_lb->end == 1)
722 pme_lb->cur = pme_lb->lower_limit;
723 pme_lb->stage = pme_lb->nstage;
725 else if (pme_lb->stage > 0 && pme_lb->end > 1)
727 /* If stage = nstage-1:
728 * scan over all setups, rerunning only those setups
729 * which are not much slower than the fastest
730 * else:
731 * use the next setup
732 * Note that we loop backward to minimize the risk of the cut-off
733 * getting limited by DD DLB, since the DLB cut-off limit is set
734 * to the fastest PME setup.
738 if (pme_lb->cur > pme_lb->start)
740 pme_lb->cur--;
742 else
744 pme_lb->stage++;
746 pme_lb->cur = pme_lb->end - 1;
749 while (pme_lb->stage == pme_lb->nstage - 1 &&
750 pme_lb->setup[pme_lb->cur].count > 0 &&
751 pme_lb->setup[pme_lb->cur].cycles > cycles_fast*maxRelativeSlowdownAccepted);
753 if (pme_lb->stage == pme_lb->nstage)
755 /* We are done optimizing, use the fastest setup we found */
756 pme_lb->cur = pme_lb->fastest;
760 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
762 OK = change_dd_cutoff(cr, state, pme_lb->setup[pme_lb->cur].rlistOuter);
763 if (!OK)
765 /* For some reason the chosen cut-off is incompatible with DD.
766 * We should continue scanning a more limited range of cut-off's.
768 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
770 /* stage=nstage says we're finished, but we should continue
771 * balancing, so we set back stage which was just incremented.
773 pme_lb->stage--;
775 if (pme_lb->cur <= pme_lb->fastest)
777 /* This should not happen, as we set limits on the DLB bounds.
778 * But we implement a complete failsafe solution anyhow.
780 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
781 "The fastest PP/PME load balancing setting (cutoff %.3d nm) is no longer available due to DD DLB or box size limitations", pme_lb->fastest);
782 pme_lb->fastest = pme_lb->lower_limit;
783 pme_lb->start = pme_lb->lower_limit;
785 /* Limit the range to below the current cut-off, scan from start */
786 pme_lb->end = pme_lb->cur;
787 pme_lb->cur = pme_lb->start;
788 pme_lb->elimited = epmelblimDD;
789 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
793 /* Change the Coulomb cut-off and the PME grid */
795 set = &pme_lb->setup[pme_lb->cur];
797 NbnxnListParameters *listParams = nbv->listParams.get();
799 ic->rcoulomb = set->rcut_coulomb;
800 listParams->rlistOuter = set->rlistOuter;
801 listParams->rlistInner = set->rlistInner;
802 ic->ewaldcoeff_q = set->ewaldcoeff_q;
803 /* TODO: centralize the code that sets the potentials shifts */
804 if (ic->coulomb_modifier == eintmodPOTSHIFT)
806 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
807 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
809 if (EVDW_PME(ic->vdwtype))
811 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
812 ic->rvdw = set->rcut_coulomb;
813 ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
814 if (ic->vdw_modifier == eintmodPOTSHIFT)
816 real crc2;
818 ic->dispersion_shift.cpot = -1.0/gmx::power6(static_cast<double>(ic->rvdw));
819 ic->repulsion_shift.cpot = -1.0/gmx::power12(static_cast<double>(ic->rvdw));
820 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
821 crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
822 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
826 /* We always re-initialize the tables whether they are used or not */
827 init_interaction_const_tables(nullptr, ic, rtab);
829 nbnxn_gpu_pme_loadbal_update_param(nbv, ic, listParams);
831 if (!pme_lb->bSepPMERanks)
833 /* FIXME:
834 * CPU PME keeps a list of allocated pmedata's, that's why pme_lb->setup[pme_lb->cur].pmedata is not always nullptr.
835 * GPU PME, however, currently needs the gmx_pme_reinit always called on load balancing
836 * (pme_gpu_reinit might be not sufficiently decoupled from gmx_pme_init).
837 * This can lead to a lot of reallocations for PME GPU.
838 * Would be nicer if the allocated grid list was hidden within a single pmedata structure.
840 if ((pme_lb->setup[pme_lb->cur].pmedata == nullptr) || pme_gpu_task_enabled(pme_lb->setup[pme_lb->cur].pmedata))
842 /* Generate a new PME data structure,
843 * copying part of the old pointers.
845 gmx_pme_reinit(&set->pmedata,
846 cr, pme_lb->setup[0].pmedata, &ir,
847 set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
849 *pmedata = set->pmedata;
851 else
853 /* Tell our PME-only rank to switch grid */
854 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
857 if (debug)
859 print_grid(nullptr, debug, "", "switched to", set, -1);
862 if (pme_lb->stage == pme_lb->nstage)
864 print_grid(fp_err, fp_log, "", "optimal", set, -1);
868 /*! \brief Prepare for another round of PME load balancing
870 * \param[in,out] pme_lb Pointer to PME load balancing struct
871 * \param[in] bDlbUnlocked TRUE is DLB was locked and is now unlocked
873 * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
874 * the PP/PME balance might change and re-balancing can improve performance.
875 * This function adds 2 stages and adjusts the considered setup range.
877 static void continue_pme_loadbal(pme_load_balancing_t *pme_lb,
878 gmx_bool bDlbUnlocked)
880 /* Add 2 tuning stages, keep the detected end of the setup range */
881 pme_lb->nstage += 2;
882 if (bDlbUnlocked && pme_lb->bSepPMERanks)
884 /* With separate PME ranks, DLB should always lower the PP load and
885 * can only increase the PME load (more communication and imbalance),
886 * so we only need to scan longer cut-off's.
888 pme_lb->lower_limit = pme_lb->cur;
890 pme_lb->start = pme_lb->lower_limit;
893 void pme_loadbal_do(pme_load_balancing_t *pme_lb,
894 t_commrec *cr,
895 FILE *fp_err,
896 FILE *fp_log,
897 const gmx::MDLogger &mdlog,
898 const t_inputrec &ir,
899 t_forcerec *fr,
900 const t_state &state,
901 gmx_wallcycle_t wcycle,
902 int64_t step,
903 int64_t step_rel,
904 gmx_bool *bPrinting)
906 int n_prev;
907 double cycles_prev;
909 assert(pme_lb != nullptr);
911 if (!pme_lb->bActive)
913 return;
916 n_prev = pme_lb->cycles_n;
917 cycles_prev = pme_lb->cycles_c;
918 wallcycle_get(wcycle, ewcSTEP, &pme_lb->cycles_n, &pme_lb->cycles_c);
920 /* Before the first step we haven't done any steps yet.
921 * Also handle cases where ir.init_step % ir.nstlist != 0.
923 if (pme_lb->cycles_n < ir.nstlist)
925 return;
927 /* Sanity check, we expect nstlist cycle counts */
928 if (pme_lb->cycles_n - n_prev != ir.nstlist)
930 /* We could return here, but it's safer to issue an error and quit */
931 gmx_incons("pme_loadbal_do called at an interval != nstlist");
934 /* PME grid + cut-off optimization with GPUs or PME ranks */
935 if (!pme_lb->bBalance && pme_lb->bSepPMERanks)
937 if (pme_lb->bTriggerOnDLB)
939 pme_lb->bBalance = dd_dlb_is_on(cr->dd);
941 /* We should ignore the first timing to avoid timing allocation
942 * overhead. And since the PME load balancing is called just
943 * before DD repartitioning, the ratio returned by dd_pme_f_ratio
944 * is not over the last nstlist steps, but the nstlist steps before
945 * that. So the first useful ratio is available at step_rel=3*nstlist.
947 else if (step_rel >= 3*ir.nstlist)
949 if (DDMASTER(cr->dd))
951 /* If PME rank load is too high, start tuning */
952 pme_lb->bBalance =
953 (dd_pme_f_ratio(cr->dd) >= loadBalanceTriggerFactor);
955 dd_bcast(cr->dd, sizeof(gmx_bool), &pme_lb->bBalance);
958 pme_lb->bActive = (pme_lb->bBalance ||
959 step_rel <= pme_lb->step_rel_stop);
962 /* The location in the code of this balancing termination is strange.
963 * You would expect to have it after the call to pme_load_balance()
964 * below, since there pme_lb->stage is updated.
965 * But when terminating directly after deciding on and selecting the
966 * optimal setup, DLB will turn on right away if it was locked before.
967 * This might be due to PME reinitialization. So we check stage here
968 * to allow for another nstlist steps with DLB locked to stabilize
969 * the performance.
971 if (pme_lb->bBalance && pme_lb->stage == pme_lb->nstage)
973 pme_lb->bBalance = FALSE;
975 if (DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
977 /* Unlock the DLB=auto, DLB is allowed to activate */
978 dd_dlb_unlock(cr->dd);
979 GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: DLB can now turn on, when beneficial");
981 /* We don't deactivate the tuning yet, since we will balance again
982 * after DLB gets turned on, if it does within PMETune_period.
984 continue_pme_loadbal(pme_lb, TRUE);
985 pme_lb->bTriggerOnDLB = TRUE;
986 pme_lb->step_rel_stop = step_rel + PMETunePeriod*ir.nstlist;
988 else
990 /* We're completely done with PME tuning */
991 pme_lb->bActive = FALSE;
994 if (DOMAINDECOMP(cr))
996 /* Set the cut-off limit to the final selected cut-off,
997 * so we don't have artificial DLB limits.
998 * This also ensures that we won't disable the currently
999 * optimal setting during a second round of PME balancing.
1001 set_dd_dlb_max_cutoff(cr, fr->nbv->listParams->rlistOuter);
1005 if (pme_lb->bBalance)
1007 /* We might not have collected nstlist steps in cycles yet,
1008 * since init_step might not be a multiple of nstlist,
1009 * but the first data collected is skipped anyhow.
1011 pme_load_balance(pme_lb, cr,
1012 fp_err, fp_log, mdlog,
1013 ir, state, pme_lb->cycles_c - cycles_prev,
1014 fr->ic, fr->nbv, &fr->pmedata,
1015 step);
1017 /* Update deprecated rlist in forcerec to stay in sync with fr->nbv */
1018 fr->rlist = fr->nbv->listParams->rlistOuter;
1020 if (ir.eDispCorr != edispcNO)
1022 calc_enervirdiff(nullptr, ir.eDispCorr, fr);
1026 if (!pme_lb->bBalance &&
1027 (!pme_lb->bSepPMERanks || step_rel > pme_lb->step_rel_stop))
1029 /* We have just deactivated the balancing and we're not measuring PP/PME
1030 * imbalance during the first steps of the run: deactivate the tuning.
1032 pme_lb->bActive = FALSE;
1035 if (!(pme_lb->bActive) && DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
1037 /* Make sure DLB is allowed when we deactivate PME tuning */
1038 dd_dlb_unlock(cr->dd);
1039 GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: DLB can now turn on, when beneficial");
1042 *bPrinting = pme_lb->bBalance;
1045 /*! \brief Return product of the number of PME grid points in each dimension */
1046 static int pme_grid_points(const pme_setup_t *setup)
1048 return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
1051 /*! \brief Print one load-balancing setting */
1052 static void print_pme_loadbal_setting(FILE *fplog,
1053 const char *name,
1054 const pme_setup_t *setup)
1056 fprintf(fplog,
1057 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
1058 name,
1059 setup->rcut_coulomb, setup->rlistInner,
1060 setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
1061 setup->spacing, 1/setup->ewaldcoeff_q);
1064 /*! \brief Print all load-balancing settings */
1065 static void print_pme_loadbal_settings(pme_load_balancing_t *pme_lb,
1066 FILE *fplog,
1067 const gmx::MDLogger &mdlog,
1068 gmx_bool bNonBondedOnGPU)
1070 double pp_ratio, grid_ratio;
1071 real pp_ratio_temporary;
1073 pp_ratio_temporary = pme_lb->setup[pme_lb->cur].rlistInner / pme_lb->setup[0].rlistInner;
1074 pp_ratio = gmx::power3(pp_ratio_temporary);
1075 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
1076 static_cast<double>(pme_grid_points(&pme_lb->setup[0]));
1078 fprintf(fplog, "\n");
1079 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
1080 fprintf(fplog, "\n");
1081 /* Here we only warn when the optimal setting is the last one */
1082 if (pme_lb->elimited != epmelblimNO &&
1083 pme_lb->cur == pme_loadbal_end(pme_lb)-1)
1085 fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
1086 pmelblim_str[pme_lb->elimited]);
1087 fprintf(fplog, " you might not have reached a good load balance.\n");
1088 if (pme_lb->elimited == epmelblimDD)
1090 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
1092 fprintf(fplog, "\n");
1094 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
1095 fprintf(fplog, " particle-particle PME\n");
1096 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
1097 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
1098 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
1099 fprintf(fplog, " cost-ratio %4.2f %4.2f\n",
1100 pp_ratio, grid_ratio);
1101 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
1103 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
1105 GMX_LOG(mdlog.warning).asParagraph().appendText(
1106 "NOTE: PME load balancing increased the non-bonded workload by more than 50%.\n"
1107 " For better performance, use (more) PME ranks (mdrun -npme),\n"
1108 " or if you are beyond the scaling limit, use fewer total ranks (or nodes).");
1110 else
1112 fprintf(fplog, "\n");
1116 void pme_loadbal_done(pme_load_balancing_t *pme_lb,
1117 FILE *fplog,
1118 const gmx::MDLogger &mdlog,
1119 gmx_bool bNonBondedOnGPU)
1121 if (fplog != nullptr && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
1123 print_pme_loadbal_settings(pme_lb, fplog, mdlog, bNonBondedOnGPU);
1126 delete pme_lb;