2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017, 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"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/fft/calcgrid.h"
60 #include "gromacs/gmxlib/md_logging.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/nbnxn_gpu_data_mgmt.h"
66 #include "gromacs/mdlib/sim_util.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/timing/wallcycle.h"
72 #include "gromacs/utility/cstringutil.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/smalloc.h"
77 #include "pme-internal.h"
79 /*! \brief Parameters and settings for one PP-PME setup */
81 real rcut_coulomb; /**< Coulomb cut-off */
82 real rlist; /**< pair-list cut-off */
83 real spacing; /**< (largest) PME grid spacing */
84 ivec grid; /**< the PME grid dimensions */
85 real grid_efficiency; /**< ineffiency factor for non-uniform grids <= 1 */
86 real ewaldcoeff_q; /**< Electrostatic Ewald coefficient */
87 real ewaldcoeff_lj; /**< LJ Ewald coefficient, only for the call to send_switchgrid */
88 struct gmx_pme_t *pmedata; /**< the data structure used in the PME code */
89 int count; /**< number of times this setup has been timed */
90 double cycles; /**< the fastest time for this setup in cycles */
93 /*! \brief After 50 nstlist periods of not observing imbalance: never tune PME */
94 const int PMETunePeriod = 50;
95 /*! \brief Trigger PME load balancing at more than 5% PME overload */
96 const real loadBalanceTriggerFactor = 1.05;
97 /*! \brief In the initial scan, step by grids that are at least a factor 0.8 coarser */
98 const real gridScaleFactor = 0.8;
99 /*! \brief In the initial scan, try to skip grids with uneven x/y/z spacing,
100 * checking if the "efficiency" is more than 5% worse than the previous grid.
102 const real relativeEfficiencyFactor = 1.05;
103 /*! \brief Rerun until a run is 12% slower setups than the fastest run so far */
104 const real maxRelativeSlowdownAccepted = 1.12;
105 /*! \brief If setups get more than 2% faster, do another round to avoid
106 * choosing a slower setup due to acceleration or fluctuations.
108 const real maxFluctuationAccepted = 1.02;
110 /*! \brief Enumeration whose values describe the effect limiting the load balancing */
112 epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimPMEGRID, epmelblimNR
115 /*! \brief Descriptive strings matching ::epmelb */
116 const char *pmelblim_str[epmelblimNR] =
117 { "no", "box size", "domain decompostion", "PME grid restriction" };
119 struct pme_load_balancing_t {
120 gmx_bool bSepPMERanks; /**< do we have separate PME ranks? */
121 gmx_bool bActive; /**< is PME tuning active? */
122 gmx_int64_t step_rel_stop; /**< stop the tuning after this value of step_rel */
123 gmx_bool bTriggerOnDLB; /**< trigger balancing only on DD DLB */
124 gmx_bool bBalance; /**< are we in the balancing phase, i.e. trying different setups? */
125 int nstage; /**< the current maximum number of stages */
127 real cut_spacing; /**< the minimum cutoff / PME grid spacing ratio */
128 real rcut_vdw; /**< Vdw cutoff (does not change) */
129 real rcut_coulomb_start; /**< Initial electrostatics cutoff */
130 real rbuf_coulomb; /**< the pairlist buffer size */
131 real rbuf_vdw; /**< the pairlist buffer size */
132 matrix box_start; /**< the initial simulation box */
133 int n; /**< the count of setup as well as the allocation size */
134 pme_setup_t *setup; /**< the PME+cutoff setups */
135 int cur; /**< the inex (in setup) of the current setup */
136 int fastest; /**< index of the fastest setup up till now */
137 int lower_limit; /**< don't go below this setup index */
138 int start; /**< start of setup index range to consider in stage>0 */
139 int end; /**< end of setup index range to consider in stage>0 */
140 int elimited; /**< was the balancing limited, uses enum above */
141 int cutoff_scheme; /**< Verlet or group cut-offs */
143 int stage; /**< the current stage */
145 int cycles_n; /**< step cycle counter cummulative count */
146 double cycles_c; /**< step cycle counter cummulative cycles */
149 /* TODO The code in this file should call this getter, rather than
150 * read bActive anywhere */
151 bool pme_loadbal_is_active(const pme_load_balancing_t *pme_lb)
153 return pme_lb != NULL && pme_lb->bActive;
156 void pme_loadbal_init(pme_load_balancing_t **pme_lb_p,
159 const t_inputrec *ir,
161 const interaction_const_t *ic,
162 struct gmx_pme_t *pmedata,
166 GMX_RELEASE_ASSERT(ir->cutoff_scheme != ecutsGROUP, "PME tuning is not supported with cutoff-scheme=group (because it contains bugs)");
168 pme_load_balancing_t *pme_lb;
172 // Note that we don't (yet) support PME load balancing with LJ-PME only.
173 GMX_RELEASE_ASSERT(EEL_PME(ir->coulombtype), "pme_loadbal_init called without PME electrostatics");
174 // To avoid complexity, we require a single cut-off with PME for q+LJ.
175 // This is checked by grompp, but it doesn't hurt to check again.
176 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");
180 pme_lb->bSepPMERanks = !(cr->duty & DUTY_PME);
182 /* Initially we turn on balancing directly on based on PP/PME imbalance */
183 pme_lb->bTriggerOnDLB = FALSE;
185 /* Any number of stages >= 2 is supported */
188 pme_lb->cutoff_scheme = ir->cutoff_scheme;
190 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
191 pme_lb->rbuf_vdw = ic->rlist - ic->rvdw;
193 copy_mat(box, pme_lb->box_start);
194 if (ir->ePBC == epbcXY && ir->nwall == 2)
196 svmul(ir->wall_ewald_zfac, pme_lb->box_start[ZZ], pme_lb->box_start[ZZ]);
200 snew(pme_lb->setup, pme_lb->n);
202 pme_lb->rcut_vdw = ic->rvdw;
203 pme_lb->rcut_coulomb_start = ir->rcoulomb;
206 pme_lb->setup[0].rcut_coulomb = ic->rcoulomb;
207 pme_lb->setup[0].rlist = ic->rlist;
208 pme_lb->setup[0].grid[XX] = ir->nkx;
209 pme_lb->setup[0].grid[YY] = ir->nky;
210 pme_lb->setup[0].grid[ZZ] = ir->nkz;
211 pme_lb->setup[0].ewaldcoeff_q = ic->ewaldcoeff_q;
212 pme_lb->setup[0].ewaldcoeff_lj = ic->ewaldcoeff_lj;
214 pme_lb->setup[0].pmedata = pmedata;
217 for (d = 0; d < DIM; d++)
219 sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
225 pme_lb->setup[0].spacing = spm;
227 if (ir->fourier_spacing > 0)
229 pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
233 pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
239 pme_lb->lower_limit = 0;
242 pme_lb->elimited = epmelblimNO;
244 pme_lb->cycles_n = 0;
245 pme_lb->cycles_c = 0;
247 if (!wallcycle_have_counter())
249 md_print_warn(cr, fp_log, "NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use PME-PP balancing.\n");
252 /* Tune with GPUs and/or separate PME ranks.
253 * When running only on a CPU without PME ranks, PME tuning will only help
254 * with small numbers of atoms in the cut-off sphere.
256 pme_lb->bActive = (wallcycle_have_counter() && (bUseGPU ||
257 pme_lb->bSepPMERanks));
259 /* With GPUs and no separate PME ranks we can't measure the PP/PME
260 * imbalance, so we start balancing right away.
261 * Otherwise we only start balancing after we observe imbalance.
263 pme_lb->bBalance = (pme_lb->bActive && (bUseGPU && !pme_lb->bSepPMERanks));
265 pme_lb->step_rel_stop = PMETunePeriod*ir->nstlist;
267 /* Delay DD load balancing when GPUs are used */
268 if (pme_lb->bActive && DOMAINDECOMP(cr) && cr->dd->nnodes > 1 && bUseGPU)
270 /* Lock DLB=auto to off (does nothing when DLB=yes/no.
271 * With GPUs and separate PME nodes, we want to first
272 * do PME tuning without DLB, since DLB might limit
273 * the cut-off, which never improves performance.
274 * We allow for DLB + PME tuning after a first round of tuning.
277 if (dd_dlb_is_locked(cr->dd))
279 md_print_warn(cr, fp_log, "NOTE: DLB will not turn on during the first phase of PME tuning\n");
285 *bPrinting = pme_lb->bBalance;
288 /*! \brief Try to increase the cutoff during load balancing */
289 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t *pme_lb,
291 const gmx_domdec_t *dd)
294 int npmeranks_x, npmeranks_y;
296 real tmpr_coulomb, tmpr_vdw;
300 /* Try to add a new setup with next larger cut-off to the list */
302 srenew(pme_lb->setup, pme_lb->n);
303 set = &pme_lb->setup[pme_lb->n-1];
306 get_pme_nnodes(dd, &npmeranks_x, &npmeranks_y);
311 /* Avoid infinite while loop, which can occur at the minimum grid size.
312 * Note that in practice load balancing will stop before this point.
313 * The factor 2.1 allows for the extreme case in which only grids
314 * of powers of 2 are allowed (the current code supports more grids).
324 clear_ivec(set->grid);
325 sp = calc_grid(NULL, pme_lb->box_start,
326 fac*pme_lb->setup[pme_lb->cur].spacing,
331 /* As here we can't easily check if one of the PME ranks
332 * uses threading, we do a conservative grid check.
333 * This means we can't use pme_order or less grid lines
334 * per PME rank along x, which is not a strong restriction.
336 gmx_pme_check_restrictions(pme_order,
337 set->grid[XX], set->grid[YY], set->grid[ZZ],
338 npmeranks_x, npmeranks_y,
343 while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
345 set->rcut_coulomb = pme_lb->cut_spacing*sp;
346 if (set->rcut_coulomb < pme_lb->rcut_coulomb_start)
348 /* This is unlikely, but can happen when e.g. continuing from
349 * a checkpoint after equilibration where the box shrank a lot.
350 * We want to avoid rcoulomb getting smaller than rvdw
351 * and there might be more issues with decreasing rcoulomb.
353 set->rcut_coulomb = pme_lb->rcut_coulomb_start;
356 if (pme_lb->cutoff_scheme == ecutsVERLET)
358 /* Never decrease the Coulomb and VdW list buffers */
359 set->rlist = std::max(set->rcut_coulomb + pme_lb->rbuf_coulomb,
360 pme_lb->rcut_vdw + pme_lb->rbuf_vdw);
364 /* TODO Remove these lines and pme_lb->cutoff_scheme */
365 tmpr_coulomb = set->rcut_coulomb + pme_lb->rbuf_coulomb;
366 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbuf_vdw;
367 /* Two (known) bugs with cutoff-scheme=group here:
368 * - This modification of rlist results in incorrect DD comunication.
369 * - We should set fr->bTwinRange = (fr->rlistlong > fr->rlist).
371 set->rlist = std::min(tmpr_coulomb, tmpr_vdw);
375 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
376 set->grid_efficiency = 1;
377 for (d = 0; d < DIM; d++)
379 set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
381 /* The Ewald coefficient is inversly proportional to the cut-off */
383 pme_lb->setup[0].ewaldcoeff_q*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
384 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
386 pme_lb->setup[0].ewaldcoeff_lj*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
393 fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
394 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
399 /*! \brief Print the PME grid */
400 static void print_grid(FILE *fp_err, FILE *fp_log,
403 const pme_setup_t *set,
406 char buf[STRLEN], buft[STRLEN];
410 sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
416 sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
418 desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
422 fprintf(fp_err, "\r%s\n", buf);
427 fprintf(fp_log, "%s\n", buf);
431 /*! \brief Return the index of the last setup used in PME load balancing */
432 static int pme_loadbal_end(pme_load_balancing_t *pme_lb)
434 /* In the initial stage only n is set; end is not set yet */
445 /*! \brief Print descriptive string about what limits PME load balancing */
446 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
448 pme_load_balancing_t *pme_lb)
450 char buf[STRLEN], sbuf[22];
452 sprintf(buf, "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
453 gmx_step_str(step, sbuf),
454 pmelblim_str[pme_lb->elimited],
455 pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
458 fprintf(fp_err, "\r%s\n", buf);
463 fprintf(fp_log, "%s\n", buf);
467 /*! \brief Switch load balancing to stage 1
469 * In this stage, only reasonably fast setups are run again. */
470 static void switch_to_stage1(pme_load_balancing_t *pme_lb)
472 /* Increase start until we find a setup that is not slower than
473 * maxRelativeSlowdownAccepted times the fastest setup.
475 pme_lb->start = pme_lb->lower_limit;
476 while (pme_lb->start + 1 < pme_lb->n &&
477 (pme_lb->setup[pme_lb->start].count == 0 ||
478 pme_lb->setup[pme_lb->start].cycles >
479 pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted))
483 /* While increasing start, we might have skipped setups that we did not
484 * time during stage 0. We want to extend the range for stage 1 to include
485 * any skipped setups that lie between setups that were measured to be
486 * acceptably fast and too slow.
488 while (pme_lb->start > pme_lb->lower_limit &&
489 pme_lb->setup[pme_lb->start - 1].count == 0)
494 /* Decrease end only with setups that we timed and that are slow. */
495 pme_lb->end = pme_lb->n;
496 if (pme_lb->setup[pme_lb->end - 1].count > 0 &&
497 pme_lb->setup[pme_lb->end - 1].cycles >
498 pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted)
505 /* Next we want to choose setup pme_lb->end-1, but as we will decrease
506 * pme_lb->cur by one right after returning, we set cur to end.
508 pme_lb->cur = pme_lb->end;
511 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
513 * The adjustment is done to generate a different non-bonded PP and PME load.
514 * With separate PME ranks (PP and PME on different processes) or with
515 * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
516 * and changing the load will affect the load balance and performance.
517 * The total time for a set of integration steps is monitored and a range
518 * of grid/cut-off setups is scanned. After calling pme_load_balance many
519 * times and acquiring enough statistics, the best performing setup is chosen.
520 * Here we try to take into account fluctuations and changes due to external
521 * factors as well as DD load balancing.
524 pme_load_balance(pme_load_balancing_t *pme_lb,
528 const t_inputrec *ir,
531 interaction_const_t *ic,
532 struct nonbonded_verlet_t *nbv,
533 struct gmx_pme_t ** pmedata,
539 char buf[STRLEN], sbuf[22];
544 gmx_sumd(1, &cycles, cr);
545 cycles /= cr->nnodes;
548 set = &pme_lb->setup[pme_lb->cur];
551 rtab = ir->rlist + ir->tabext;
553 if (set->count % 2 == 1)
555 /* Skip the first cycle, because the first step after a switch
556 * is much slower due to allocation and/or caching effects.
561 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
562 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
566 set->cycles = cycles;
570 if (cycles*maxFluctuationAccepted < set->cycles &&
571 pme_lb->stage == pme_lb->nstage - 1)
573 /* The performance went up a lot (due to e.g. DD load balancing).
574 * Add a stage, keep the minima, but rescan all setups.
580 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
581 "Increased the number stages to %d"
582 " and ignoring the previous performance\n",
583 set->grid[XX], set->grid[YY], set->grid[ZZ],
584 set->cycles*1e-6, cycles*1e-6, maxFluctuationAccepted,
588 set->cycles = std::min(set->cycles, cycles);
591 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
593 pme_lb->fastest = pme_lb->cur;
595 if (DOMAINDECOMP(cr))
597 /* We found a new fastest setting, ensure that with subsequent
598 * shorter cut-off's the dynamic load balancing does not make
599 * the use of the current cut-off impossible. This solution is
600 * a trade-off, as the PME load balancing and DD domain size
601 * load balancing can interact in complex ways.
602 * With the Verlet kernels, DD load imbalance will usually be
603 * mainly due to bonded interaction imbalance, which will often
604 * quickly push the domain boundaries beyond the limit for the
605 * optimal, PME load balanced, cut-off. But it could be that
606 * better overal performance can be obtained with a slightly
607 * shorter cut-off and better DD load balancing.
609 set_dd_dlb_max_cutoff(cr, pme_lb->setup[pme_lb->fastest].rlist);
612 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
614 /* Check in stage 0 if we should stop scanning grids.
615 * Stop when the time is more than maxRelativeSlowDownAccepted longer than the fastest.
617 if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
618 cycles > pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted)
620 pme_lb->n = pme_lb->cur + 1;
621 /* Done with scanning, go to stage 1 */
622 switch_to_stage1(pme_lb);
625 if (pme_lb->stage == 0)
629 gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
633 if (pme_lb->cur+1 < pme_lb->n)
635 /* We had already generated the next setup */
640 /* Find the next setup */
641 OK = pme_loadbal_increase_cutoff(pme_lb, ir->pme_order, cr->dd);
645 pme_lb->elimited = epmelblimPMEGRID;
649 if (OK && ir->ePBC != epbcNONE)
651 OK = (gmx::square(pme_lb->setup[pme_lb->cur+1].rlist)
652 <= max_cutoff2(ir->ePBC, state->box));
655 pme_lb->elimited = epmelblimBOX;
663 if (DOMAINDECOMP(cr))
665 OK = change_dd_cutoff(cr, state, ir,
666 pme_lb->setup[pme_lb->cur].rlist);
669 /* Failed: do not use this setup */
671 pme_lb->elimited = epmelblimDD;
677 /* We hit the upper limit for the cut-off,
678 * the setup should not go further than cur.
680 pme_lb->n = pme_lb->cur + 1;
681 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
682 /* Switch to the next stage */
683 switch_to_stage1(pme_lb);
687 !(pme_lb->setup[pme_lb->cur].grid[XX]*
688 pme_lb->setup[pme_lb->cur].grid[YY]*
689 pme_lb->setup[pme_lb->cur].grid[ZZ] <
690 gridsize_start*gridScaleFactor
692 pme_lb->setup[pme_lb->cur].grid_efficiency <
693 pme_lb->setup[pme_lb->cur-1].grid_efficiency*relativeEfficiencyFactor));
696 if (pme_lb->stage > 0 && pme_lb->end == 1)
698 pme_lb->cur = pme_lb->lower_limit;
699 pme_lb->stage = pme_lb->nstage;
701 else if (pme_lb->stage > 0 && pme_lb->end > 1)
703 /* If stage = nstage-1:
704 * scan over all setups, rerunning only those setups
705 * which are not much slower than the fastest
708 * Note that we loop backward to minimize the risk of the cut-off
709 * getting limited by DD DLB, since the DLB cut-off limit is set
710 * to the fastest PME setup.
714 if (pme_lb->cur > pme_lb->start)
722 pme_lb->cur = pme_lb->end - 1;
725 while (pme_lb->stage == pme_lb->nstage - 1 &&
726 pme_lb->setup[pme_lb->cur].count > 0 &&
727 pme_lb->setup[pme_lb->cur].cycles > cycles_fast*maxRelativeSlowdownAccepted);
729 if (pme_lb->stage == pme_lb->nstage)
731 /* We are done optimizing, use the fastest setup we found */
732 pme_lb->cur = pme_lb->fastest;
736 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
738 OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlist);
741 /* For some reason the chosen cut-off is incompatible with DD.
742 * We should continue scanning a more limited range of cut-off's.
744 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
746 /* stage=nstage says we're finished, but we should continue
747 * balancing, so we set back stage which was just incremented.
751 if (pme_lb->cur <= pme_lb->fastest)
753 /* This should not happen, as we set limits on the DLB bounds.
754 * But we implement a complete failsafe solution anyhow.
756 md_print_warn(cr, fp_log, "The fastest PP/PME load balancing setting (cutoff %.3f nm) is no longer available due to DD DLB or box size limitations\n");
757 pme_lb->fastest = pme_lb->lower_limit;
758 pme_lb->start = pme_lb->lower_limit;
760 /* Limit the range to below the current cut-off, scan from start */
761 pme_lb->end = pme_lb->cur;
762 pme_lb->cur = pme_lb->start;
763 pme_lb->elimited = epmelblimDD;
764 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
768 /* Change the Coulomb cut-off and the PME grid */
770 set = &pme_lb->setup[pme_lb->cur];
772 ic->rcoulomb = set->rcut_coulomb;
773 ic->rlist = set->rlist;
774 ic->ewaldcoeff_q = set->ewaldcoeff_q;
775 /* TODO: centralize the code that sets the potentials shifts */
776 if (ic->coulomb_modifier == eintmodPOTSHIFT)
778 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
779 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
781 if (EVDW_PME(ic->vdwtype))
783 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
784 ic->rvdw = set->rcut_coulomb;
785 ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
786 if (ic->vdw_modifier == eintmodPOTSHIFT)
790 ic->dispersion_shift.cpot = -1.0/gmx::power6(static_cast<double>(ic->rvdw));
791 ic->repulsion_shift.cpot = -1.0/gmx::power12(static_cast<double>(ic->rvdw));
792 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
793 crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
794 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
798 /* We always re-initialize the tables whether they are used or not */
799 init_interaction_const_tables(NULL, ic, rtab);
801 nbnxn_gpu_pme_loadbal_update_param(nbv, ic);
803 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
804 * also sharing texture references. To keep the code simple, we don't
805 * treat texture references as shared resources, but this means that
806 * the coulomb_tab texture ref will get updated by multiple threads.
807 * Hence, to ensure that the non-bonded kernels don't start before all
808 * texture binding operations are finished, we need to wait for all ranks
809 * to arrive here before continuing.
811 * Note that we could omit this barrier if GPUs are not shared (or
812 * texture objects are used), but as this is initialization code, there
813 * is not point in complicating things.
816 if (PAR(cr) && use_GPU(nbv))
820 #endif /* GMX_THREAD_MPI */
822 if (!pme_lb->bSepPMERanks)
824 if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
826 /* Generate a new PME data structure,
827 * copying part of the old pointers.
829 gmx_pme_reinit(&set->pmedata,
830 cr, pme_lb->setup[0].pmedata, ir,
833 *pmedata = set->pmedata;
837 /* Tell our PME-only rank to switch grid */
838 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
843 print_grid(NULL, debug, "", "switched to", set, -1);
846 if (pme_lb->stage == pme_lb->nstage)
848 print_grid(fp_err, fp_log, "", "optimal", set, -1);
852 /*! \brief Prepare for another round of PME load balancing
854 * \param[in,out] pme_lb Pointer to PME load balancing struct
855 * \param[in] bDlbUnlocked TRUE is DLB was locked and is now unlocked
857 * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
858 * the PP/PME balance might change and re-balancing can improve performance.
859 * This function adds 2 stages and adjusts the considered setup range.
861 static void continue_pme_loadbal(pme_load_balancing_t *pme_lb,
862 gmx_bool bDlbUnlocked)
864 /* Add 2 tuning stages, keep the detected end of the setup range */
866 if (bDlbUnlocked && pme_lb->bSepPMERanks)
868 /* With separate PME ranks, DLB should always lower the PP load and
869 * can only increase the PME load (more communication and imbalance),
870 * so we only need to scan longer cut-off's.
872 pme_lb->lower_limit = pme_lb->cur;
874 pme_lb->start = pme_lb->lower_limit;
877 void pme_loadbal_do(pme_load_balancing_t *pme_lb,
881 const t_inputrec *ir,
884 gmx_wallcycle_t wcycle,
886 gmx_int64_t step_rel,
892 assert(pme_lb != NULL);
894 if (!pme_lb->bActive)
899 n_prev = pme_lb->cycles_n;
900 cycles_prev = pme_lb->cycles_c;
901 wallcycle_get(wcycle, ewcSTEP, &pme_lb->cycles_n, &pme_lb->cycles_c);
903 /* Before the first step we haven't done any steps yet.
904 * Also handle cases where ir->init_step % ir->nstlist != 0.
906 if (pme_lb->cycles_n < ir->nstlist)
910 /* Sanity check, we expect nstlist cycle counts */
911 if (pme_lb->cycles_n - n_prev != ir->nstlist)
913 /* We could return here, but it's safer to issue an error and quit */
914 gmx_incons("pme_loadbal_do called at an interval != nstlist");
917 /* PME grid + cut-off optimization with GPUs or PME ranks */
918 if (!pme_lb->bBalance && pme_lb->bSepPMERanks)
920 if (pme_lb->bTriggerOnDLB)
922 pme_lb->bBalance = dd_dlb_is_on(cr->dd);
924 /* We should ignore the first timing to avoid timing allocation
925 * overhead. And since the PME load balancing is called just
926 * before DD repartitioning, the ratio returned by dd_pme_f_ratio
927 * is not over the last nstlist steps, but the nstlist steps before
928 * that. So the first useful ratio is available at step_rel=3*nstlist.
930 else if (step_rel >= 3*ir->nstlist)
932 if (DDMASTER(cr->dd))
934 /* If PME rank load is too high, start tuning */
936 (dd_pme_f_ratio(cr->dd) >= loadBalanceTriggerFactor);
938 dd_bcast(cr->dd, sizeof(gmx_bool), &pme_lb->bBalance);
941 pme_lb->bActive = (pme_lb->bBalance ||
942 step_rel <= pme_lb->step_rel_stop);
945 /* The location in the code of this balancing termination is strange.
946 * You would expect to have it after the call to pme_load_balance()
947 * below, since there pme_lb->stage is updated.
948 * But when terminating directly after deciding on and selecting the
949 * optimal setup, DLB will turn on right away if it was locked before.
950 * This might be due to PME reinitialization. So we check stage here
951 * to allow for another nstlist steps with DLB locked to stabilize
954 if (pme_lb->bBalance && pme_lb->stage == pme_lb->nstage)
956 pme_lb->bBalance = FALSE;
958 if (DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
960 /* Unlock the DLB=auto, DLB is allowed to activate */
961 dd_dlb_unlock(cr->dd);
962 md_print_warn(cr, fp_log, "NOTE: DLB can now turn on, when beneficial\n");
964 /* We don't deactivate the tuning yet, since we will balance again
965 * after DLB gets turned on, if it does within PMETune_period.
967 continue_pme_loadbal(pme_lb, TRUE);
968 pme_lb->bTriggerOnDLB = TRUE;
969 pme_lb->step_rel_stop = step_rel + PMETunePeriod*ir->nstlist;
973 /* We're completely done with PME tuning */
974 pme_lb->bActive = FALSE;
977 if (DOMAINDECOMP(cr))
979 /* Set the cut-off limit to the final selected cut-off,
980 * so we don't have artificial DLB limits.
981 * This also ensures that we won't disable the currently
982 * optimal setting during a second round of PME balancing.
984 set_dd_dlb_max_cutoff(cr, fr->ic->rlist);
988 if (pme_lb->bBalance)
990 /* We might not have collected nstlist steps in cycles yet,
991 * since init_step might not be a multiple of nstlist,
992 * but the first data collected is skipped anyhow.
994 pme_load_balance(pme_lb, cr,
996 ir, state, pme_lb->cycles_c - cycles_prev,
997 fr->ic, fr->nbv, &fr->pmedata,
1000 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1001 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1002 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1003 fr->rlist = fr->ic->rlist;
1004 fr->rcoulomb = fr->ic->rcoulomb;
1005 fr->rvdw = fr->ic->rvdw;
1007 if (ir->eDispCorr != edispcNO)
1009 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1013 if (!pme_lb->bBalance &&
1014 (!pme_lb->bSepPMERanks || step_rel > pme_lb->step_rel_stop))
1016 /* We have just deactivated the balancing and we're not measuring PP/PME
1017 * imbalance during the first steps of the run: deactivate the tuning.
1019 pme_lb->bActive = FALSE;
1022 if (!(pme_lb->bActive) && DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
1024 /* Make sure DLB is allowed when we deactivate PME tuning */
1025 dd_dlb_unlock(cr->dd);
1026 md_print_warn(cr, fp_log, "NOTE: DLB can now turn on, when beneficial\n");
1029 *bPrinting = pme_lb->bBalance;
1032 /*! \brief Return product of the number of PME grid points in each dimension */
1033 static int pme_grid_points(const pme_setup_t *setup)
1035 return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
1038 /*! \brief Print one load-balancing setting */
1039 static void print_pme_loadbal_setting(FILE *fplog,
1041 const pme_setup_t *setup)
1044 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
1046 setup->rcut_coulomb, setup->rlist,
1047 setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
1048 setup->spacing, 1/setup->ewaldcoeff_q);
1051 /*! \brief Print all load-balancing settings */
1052 static void print_pme_loadbal_settings(pme_load_balancing_t *pme_lb,
1055 gmx_bool bNonBondedOnGPU)
1057 double pp_ratio, grid_ratio;
1058 real pp_ratio_temporary;
1060 pp_ratio_temporary = pme_lb->setup[pme_lb->cur].rlist / pme_lb->setup[0].rlist;
1061 pp_ratio = gmx::power3(pp_ratio_temporary);
1062 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
1063 (double)pme_grid_points(&pme_lb->setup[0]);
1065 fprintf(fplog, "\n");
1066 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
1067 fprintf(fplog, "\n");
1068 /* Here we only warn when the optimal setting is the last one */
1069 if (pme_lb->elimited != epmelblimNO &&
1070 pme_lb->cur == pme_loadbal_end(pme_lb)-1)
1072 fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
1073 pmelblim_str[pme_lb->elimited]);
1074 fprintf(fplog, " you might not have reached a good load balance.\n");
1075 if (pme_lb->elimited == epmelblimDD)
1077 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
1079 fprintf(fplog, "\n");
1081 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
1082 fprintf(fplog, " particle-particle PME\n");
1083 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
1084 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
1085 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
1086 fprintf(fplog, " cost-ratio %4.2f %4.2f\n",
1087 pp_ratio, grid_ratio);
1088 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
1090 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
1092 md_print_warn(cr, fplog,
1093 "NOTE: PME load balancing increased the non-bonded workload by more than 50%%.\n"
1094 " For better performance, use (more) PME ranks (mdrun -npme),\n"
1095 " or if you are beyond the scaling limit, use fewer total ranks (or nodes).\n");
1099 fprintf(fplog, "\n");
1103 void pme_loadbal_done(pme_load_balancing_t *pme_lb,
1106 gmx_bool bNonBondedOnGPU)
1108 if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
1110 print_pme_loadbal_settings(pme_lb, cr, fplog, bNonBondedOnGPU);
1113 /* TODO: Here we should free all pointers in pme_lb,
1114 * but as it contains pme data structures,
1115 * we need to first make pme.c free all data.