2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016, 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 pme_load_balancing_t *pme_lb;
170 // Note that we don't (yet) support PME load balancing with LJ-PME only.
171 GMX_RELEASE_ASSERT(EEL_PME(ir->coulombtype), "pme_loadbal_init called without PME electrostatics");
172 // To avoid complexity, we require a single cut-off with PME for q+LJ.
173 // This is checked by grompp, but it doesn't hurt to check again.
174 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");
178 pme_lb->bSepPMERanks = !(cr->duty & DUTY_PME);
180 /* Initially we turn on balancing directly on based on PP/PME imbalance */
181 pme_lb->bTriggerOnDLB = FALSE;
183 /* Any number of stages >= 2 is supported */
186 pme_lb->cutoff_scheme = ir->cutoff_scheme;
188 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
189 pme_lb->rbuf_vdw = ic->rlist - ic->rvdw;
191 copy_mat(box, pme_lb->box_start);
192 if (ir->ePBC == epbcXY && ir->nwall == 2)
194 svmul(ir->wall_ewald_zfac, pme_lb->box_start[ZZ], pme_lb->box_start[ZZ]);
198 snew(pme_lb->setup, pme_lb->n);
200 pme_lb->rcut_vdw = ic->rvdw;
201 pme_lb->rcut_coulomb_start = ir->rcoulomb;
204 pme_lb->setup[0].rcut_coulomb = ic->rcoulomb;
205 pme_lb->setup[0].rlist = ic->rlist;
206 pme_lb->setup[0].grid[XX] = ir->nkx;
207 pme_lb->setup[0].grid[YY] = ir->nky;
208 pme_lb->setup[0].grid[ZZ] = ir->nkz;
209 pme_lb->setup[0].ewaldcoeff_q = ic->ewaldcoeff_q;
210 pme_lb->setup[0].ewaldcoeff_lj = ic->ewaldcoeff_lj;
212 pme_lb->setup[0].pmedata = pmedata;
215 for (d = 0; d < DIM; d++)
217 sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
223 pme_lb->setup[0].spacing = spm;
225 if (ir->fourier_spacing > 0)
227 pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
231 pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
237 pme_lb->lower_limit = 0;
240 pme_lb->elimited = epmelblimNO;
242 pme_lb->cycles_n = 0;
243 pme_lb->cycles_c = 0;
245 if (!wallcycle_have_counter())
247 md_print_warn(cr, fp_log, "NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use PME-PP balancing.\n");
250 /* Tune with GPUs and/or separate PME ranks.
251 * When running only on a CPU without PME ranks, PME tuning will only help
252 * with small numbers of atoms in the cut-off sphere.
254 pme_lb->bActive = (wallcycle_have_counter() && (bUseGPU ||
255 pme_lb->bSepPMERanks));
257 /* With GPUs and no separate PME ranks we can't measure the PP/PME
258 * imbalance, so we start balancing right away.
259 * Otherwise we only start balancing after we observe imbalance.
261 pme_lb->bBalance = (pme_lb->bActive && (bUseGPU && !pme_lb->bSepPMERanks));
263 pme_lb->step_rel_stop = PMETunePeriod*ir->nstlist;
265 /* Delay DD load balancing when GPUs are used */
266 if (pme_lb->bActive && DOMAINDECOMP(cr) && cr->dd->nnodes > 1 && bUseGPU)
268 /* Lock DLB=auto to off (does nothing when DLB=yes/no.
269 * With GPUs and separate PME nodes, we want to first
270 * do PME tuning without DLB, since DLB might limit
271 * the cut-off, which never improves performance.
272 * We allow for DLB + PME tuning after a first round of tuning.
275 if (dd_dlb_is_locked(cr->dd))
277 md_print_warn(cr, fp_log, "NOTE: DLB will not turn on during the first phase of PME tuning\n");
283 *bPrinting = pme_lb->bBalance;
286 /*! \brief Try to increase the cutoff during load balancing */
287 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t *pme_lb,
289 const gmx_domdec_t *dd)
292 int npmeranks_x, npmeranks_y;
294 real tmpr_coulomb, tmpr_vdw;
298 /* Try to add a new setup with next larger cut-off to the list */
300 srenew(pme_lb->setup, pme_lb->n);
301 set = &pme_lb->setup[pme_lb->n-1];
304 get_pme_nnodes(dd, &npmeranks_x, &npmeranks_y);
309 /* Avoid infinite while loop, which can occur at the minimum grid size.
310 * Note that in practice load balancing will stop before this point.
311 * The factor 2.1 allows for the extreme case in which only grids
312 * of powers of 2 are allowed (the current code supports more grids).
322 clear_ivec(set->grid);
323 sp = calc_grid(NULL, pme_lb->box_start,
324 fac*pme_lb->setup[pme_lb->cur].spacing,
329 /* As here we can't easily check if one of the PME ranks
330 * uses threading, we do a conservative grid check.
331 * This means we can't use pme_order or less grid lines
332 * per PME rank along x, which is not a strong restriction.
334 gmx_pme_check_restrictions(pme_order,
335 set->grid[XX], set->grid[YY], set->grid[ZZ],
336 npmeranks_x, npmeranks_y,
341 while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
343 set->rcut_coulomb = pme_lb->cut_spacing*sp;
344 if (set->rcut_coulomb < pme_lb->rcut_coulomb_start)
346 /* This is unlikely, but can happen when e.g. continuing from
347 * a checkpoint after equilibration where the box shrank a lot.
348 * We want to avoid rcoulomb getting smaller than rvdw
349 * and there might be more issues with decreasing rcoulomb.
351 set->rcut_coulomb = pme_lb->rcut_coulomb_start;
354 if (pme_lb->cutoff_scheme == ecutsVERLET)
356 /* Never decrease the Coulomb and VdW list buffers */
357 set->rlist = std::max(set->rcut_coulomb + pme_lb->rbuf_coulomb,
358 pme_lb->rcut_vdw + pme_lb->rbuf_vdw);
362 tmpr_coulomb = set->rcut_coulomb + pme_lb->rbuf_coulomb;
363 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbuf_vdw;
364 set->rlist = std::min(tmpr_coulomb, tmpr_vdw);
368 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
369 set->grid_efficiency = 1;
370 for (d = 0; d < DIM; d++)
372 set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
374 /* The Ewald coefficient is inversly proportional to the cut-off */
376 pme_lb->setup[0].ewaldcoeff_q*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
377 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
379 pme_lb->setup[0].ewaldcoeff_lj*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
386 fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
387 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
392 /*! \brief Print the PME grid */
393 static void print_grid(FILE *fp_err, FILE *fp_log,
396 const pme_setup_t *set,
399 char buf[STRLEN], buft[STRLEN];
403 sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
409 sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
411 desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
415 fprintf(fp_err, "\r%s\n", buf);
420 fprintf(fp_log, "%s\n", buf);
424 /*! \brief Return the index of the last setup used in PME load balancing */
425 static int pme_loadbal_end(pme_load_balancing_t *pme_lb)
427 /* In the initial stage only n is set; end is not set yet */
438 /*! \brief Print descriptive string about what limits PME load balancing */
439 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
441 pme_load_balancing_t *pme_lb)
443 char buf[STRLEN], sbuf[22];
445 sprintf(buf, "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
446 gmx_step_str(step, sbuf),
447 pmelblim_str[pme_lb->elimited],
448 pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
451 fprintf(fp_err, "\r%s\n", buf);
456 fprintf(fp_log, "%s\n", buf);
460 /*! \brief Switch load balancing to stage 1
462 * In this stage, only reasonably fast setups are run again. */
463 static void switch_to_stage1(pme_load_balancing_t *pme_lb)
465 /* Increase start until we find a setup that is not slower than
466 * maxRelativeSlowdownAccepted times the fastest setup.
468 pme_lb->start = pme_lb->lower_limit;
469 while (pme_lb->start + 1 < pme_lb->n &&
470 (pme_lb->setup[pme_lb->start].count == 0 ||
471 pme_lb->setup[pme_lb->start].cycles >
472 pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted))
476 /* While increasing start, we might have skipped setups that we did not
477 * time during stage 0. We want to extend the range for stage 1 to include
478 * any skipped setups that lie between setups that were measured to be
479 * acceptably fast and too slow.
481 while (pme_lb->start > pme_lb->lower_limit &&
482 pme_lb->setup[pme_lb->start - 1].count == 0)
487 /* Decrease end only with setups that we timed and that are slow. */
488 pme_lb->end = pme_lb->n;
489 if (pme_lb->setup[pme_lb->end - 1].count > 0 &&
490 pme_lb->setup[pme_lb->end - 1].cycles >
491 pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted)
498 /* Next we want to choose setup pme_lb->end-1, but as we will decrease
499 * pme_lb->cur by one right after returning, we set cur to end.
501 pme_lb->cur = pme_lb->end;
504 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
506 * The adjustment is done to generate a different non-bonded PP and PME load.
507 * With separate PME ranks (PP and PME on different processes) or with
508 * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
509 * and changing the load will affect the load balance and performance.
510 * The total time for a set of integration steps is monitored and a range
511 * of grid/cut-off setups is scanned. After calling pme_load_balance many
512 * times and acquiring enough statistics, the best performing setup is chosen.
513 * Here we try to take into account fluctuations and changes due to external
514 * factors as well as DD load balancing.
517 pme_load_balance(pme_load_balancing_t *pme_lb,
521 const t_inputrec *ir,
524 interaction_const_t *ic,
525 struct nonbonded_verlet_t *nbv,
526 struct gmx_pme_t ** pmedata,
532 char buf[STRLEN], sbuf[22];
537 gmx_sumd(1, &cycles, cr);
538 cycles /= cr->nnodes;
541 set = &pme_lb->setup[pme_lb->cur];
544 rtab = ir->rlist + ir->tabext;
546 if (set->count % 2 == 1)
548 /* Skip the first cycle, because the first step after a switch
549 * is much slower due to allocation and/or caching effects.
554 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
555 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
559 set->cycles = cycles;
563 if (cycles*maxFluctuationAccepted < set->cycles &&
564 pme_lb->stage == pme_lb->nstage - 1)
566 /* The performance went up a lot (due to e.g. DD load balancing).
567 * Add a stage, keep the minima, but rescan all setups.
573 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
574 "Increased the number stages to %d"
575 " and ignoring the previous performance\n",
576 set->grid[XX], set->grid[YY], set->grid[ZZ],
577 set->cycles*1e-6, cycles*1e-6, maxFluctuationAccepted,
581 set->cycles = std::min(set->cycles, cycles);
584 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
586 pme_lb->fastest = pme_lb->cur;
588 if (DOMAINDECOMP(cr))
590 /* We found a new fastest setting, ensure that with subsequent
591 * shorter cut-off's the dynamic load balancing does not make
592 * the use of the current cut-off impossible. This solution is
593 * a trade-off, as the PME load balancing and DD domain size
594 * load balancing can interact in complex ways.
595 * With the Verlet kernels, DD load imbalance will usually be
596 * mainly due to bonded interaction imbalance, which will often
597 * quickly push the domain boundaries beyond the limit for the
598 * optimal, PME load balanced, cut-off. But it could be that
599 * better overal performance can be obtained with a slightly
600 * shorter cut-off and better DD load balancing.
602 set_dd_dlb_max_cutoff(cr, pme_lb->setup[pme_lb->fastest].rlist);
605 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
607 /* Check in stage 0 if we should stop scanning grids.
608 * Stop when the time is more than maxRelativeSlowDownAccepted longer than the fastest.
610 if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
611 cycles > pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted)
613 pme_lb->n = pme_lb->cur + 1;
614 /* Done with scanning, go to stage 1 */
615 switch_to_stage1(pme_lb);
618 if (pme_lb->stage == 0)
622 gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
626 if (pme_lb->cur+1 < pme_lb->n)
628 /* We had already generated the next setup */
633 /* Find the next setup */
634 OK = pme_loadbal_increase_cutoff(pme_lb, ir->pme_order, cr->dd);
638 pme_lb->elimited = epmelblimPMEGRID;
642 if (OK && ir->ePBC != epbcNONE)
644 OK = (gmx::square(pme_lb->setup[pme_lb->cur+1].rlist)
645 <= max_cutoff2(ir->ePBC, state->box));
648 pme_lb->elimited = epmelblimBOX;
656 if (DOMAINDECOMP(cr))
658 OK = change_dd_cutoff(cr, state, ir,
659 pme_lb->setup[pme_lb->cur].rlist);
662 /* Failed: do not use this setup */
664 pme_lb->elimited = epmelblimDD;
670 /* We hit the upper limit for the cut-off,
671 * the setup should not go further than cur.
673 pme_lb->n = pme_lb->cur + 1;
674 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
675 /* Switch to the next stage */
676 switch_to_stage1(pme_lb);
680 !(pme_lb->setup[pme_lb->cur].grid[XX]*
681 pme_lb->setup[pme_lb->cur].grid[YY]*
682 pme_lb->setup[pme_lb->cur].grid[ZZ] <
683 gridsize_start*gridScaleFactor
685 pme_lb->setup[pme_lb->cur].grid_efficiency <
686 pme_lb->setup[pme_lb->cur-1].grid_efficiency*relativeEfficiencyFactor));
689 if (pme_lb->stage > 0 && pme_lb->end == 1)
691 pme_lb->cur = pme_lb->lower_limit;
692 pme_lb->stage = pme_lb->nstage;
694 else if (pme_lb->stage > 0 && pme_lb->end > 1)
696 /* If stage = nstage-1:
697 * scan over all setups, rerunning only those setups
698 * which are not much slower than the fastest
701 * Note that we loop backward to minimize the risk of the cut-off
702 * getting limited by DD DLB, since the DLB cut-off limit is set
703 * to the fastest PME setup.
707 if (pme_lb->cur > pme_lb->start)
715 pme_lb->cur = pme_lb->end - 1;
718 while (pme_lb->stage == pme_lb->nstage - 1 &&
719 pme_lb->setup[pme_lb->cur].count > 0 &&
720 pme_lb->setup[pme_lb->cur].cycles > cycles_fast*maxRelativeSlowdownAccepted);
722 if (pme_lb->stage == pme_lb->nstage)
724 /* We are done optimizing, use the fastest setup we found */
725 pme_lb->cur = pme_lb->fastest;
729 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
731 OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlist);
734 /* For some reason the chosen cut-off is incompatible with DD.
735 * We should continue scanning a more limited range of cut-off's.
737 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
739 /* stage=nstage says we're finished, but we should continue
740 * balancing, so we set back stage which was just incremented.
744 if (pme_lb->cur <= pme_lb->fastest)
746 /* This should not happen, as we set limits on the DLB bounds.
747 * But we implement a complete failsafe solution anyhow.
749 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");
750 pme_lb->fastest = pme_lb->lower_limit;
751 pme_lb->start = pme_lb->lower_limit;
753 /* Limit the range to below the current cut-off, scan from start */
754 pme_lb->end = pme_lb->cur;
755 pme_lb->cur = pme_lb->start;
756 pme_lb->elimited = epmelblimDD;
757 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
761 /* Change the Coulomb cut-off and the PME grid */
763 set = &pme_lb->setup[pme_lb->cur];
765 ic->rcoulomb = set->rcut_coulomb;
766 ic->rlist = set->rlist;
767 ic->ewaldcoeff_q = set->ewaldcoeff_q;
768 /* TODO: centralize the code that sets the potentials shifts */
769 if (ic->coulomb_modifier == eintmodPOTSHIFT)
771 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb);
773 if (EVDW_PME(ic->vdwtype))
775 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
776 ic->rvdw = set->rcut_coulomb;
777 ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
778 if (ic->vdw_modifier == eintmodPOTSHIFT)
782 ic->dispersion_shift.cpot = -1.0/gmx::power6(static_cast<double>(ic->rvdw));
783 ic->repulsion_shift.cpot = -1.0/gmx::power12(static_cast<double>(ic->rvdw));
784 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
785 crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
786 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
790 /* We always re-initialize the tables whether they are used or not */
791 init_interaction_const_tables(NULL, ic, rtab);
793 nbnxn_gpu_pme_loadbal_update_param(nbv, ic);
795 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
796 * also sharing texture references. To keep the code simple, we don't
797 * treat texture references as shared resources, but this means that
798 * the coulomb_tab texture ref will get updated by multiple threads.
799 * Hence, to ensure that the non-bonded kernels don't start before all
800 * texture binding operations are finished, we need to wait for all ranks
801 * to arrive here before continuing.
803 * Note that we could omit this barrier if GPUs are not shared (or
804 * texture objects are used), but as this is initialization code, there
805 * is not point in complicating things.
808 if (PAR(cr) && use_GPU(nbv))
812 #endif /* GMX_THREAD_MPI */
814 if (!pme_lb->bSepPMERanks)
816 if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
818 /* Generate a new PME data structure,
819 * copying part of the old pointers.
821 gmx_pme_reinit(&set->pmedata,
822 cr, pme_lb->setup[0].pmedata, ir,
825 *pmedata = set->pmedata;
829 /* Tell our PME-only rank to switch grid */
830 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
835 print_grid(NULL, debug, "", "switched to", set, -1);
838 if (pme_lb->stage == pme_lb->nstage)
840 print_grid(fp_err, fp_log, "", "optimal", set, -1);
844 /*! \brief Prepare for another round of PME load balancing
846 * \param[in,out] pme_lb Pointer to PME load balancing struct
847 * \param[in] bDlbUnlocked TRUE is DLB was locked and is now unlocked
849 * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
850 * the PP/PME balance might change and re-balancing can improve performance.
851 * This function adds 2 stages and adjusts the considered setup range.
853 static void continue_pme_loadbal(pme_load_balancing_t *pme_lb,
854 gmx_bool bDlbUnlocked)
856 /* Add 2 tuning stages, keep the detected end of the setup range */
858 if (bDlbUnlocked && pme_lb->bSepPMERanks)
860 /* With separate PME ranks, DLB should always lower the PP load and
861 * can only increase the PME load (more communication and imbalance),
862 * so we only need to scan longer cut-off's.
864 pme_lb->lower_limit = pme_lb->cur;
866 pme_lb->start = pme_lb->lower_limit;
869 void pme_loadbal_do(pme_load_balancing_t *pme_lb,
873 const t_inputrec *ir,
876 gmx_wallcycle_t wcycle,
878 gmx_int64_t step_rel,
884 assert(pme_lb != NULL);
886 if (!pme_lb->bActive)
891 n_prev = pme_lb->cycles_n;
892 cycles_prev = pme_lb->cycles_c;
893 wallcycle_get(wcycle, ewcSTEP, &pme_lb->cycles_n, &pme_lb->cycles_c);
895 /* Before the first step we haven't done any steps yet.
896 * Also handle cases where ir->init_step % ir->nstlist != 0.
898 if (pme_lb->cycles_n < ir->nstlist)
902 /* Sanity check, we expect nstlist cycle counts */
903 if (pme_lb->cycles_n - n_prev != ir->nstlist)
905 /* We could return here, but it's safer to issue an error and quit */
906 gmx_incons("pme_loadbal_do called at an interval != nstlist");
909 /* PME grid + cut-off optimization with GPUs or PME ranks */
910 if (!pme_lb->bBalance && pme_lb->bSepPMERanks)
912 if (pme_lb->bTriggerOnDLB)
914 pme_lb->bBalance = dd_dlb_is_on(cr->dd);
916 /* We should ignore the first timing to avoid timing allocation
917 * overhead. And since the PME load balancing is called just
918 * before DD repartitioning, the ratio returned by dd_pme_f_ratio
919 * is not over the last nstlist steps, but the nstlist steps before
920 * that. So the first useful ratio is available at step_rel=3*nstlist.
922 else if (step_rel >= 3*ir->nstlist)
924 if (DDMASTER(cr->dd))
926 /* If PME rank load is too high, start tuning */
928 (dd_pme_f_ratio(cr->dd) >= loadBalanceTriggerFactor);
930 dd_bcast(cr->dd, sizeof(gmx_bool), &pme_lb->bBalance);
933 pme_lb->bActive = (pme_lb->bBalance ||
934 step_rel <= pme_lb->step_rel_stop);
937 /* The location in the code of this balancing termination is strange.
938 * You would expect to have it after the call to pme_load_balance()
939 * below, since there pme_lb->stage is updated.
940 * But when terminating directly after deciding on and selecting the
941 * optimal setup, DLB will turn on right away if it was locked before.
942 * This might be due to PME reinitialization. So we check stage here
943 * to allow for another nstlist steps with DLB locked to stabilize
946 if (pme_lb->bBalance && pme_lb->stage == pme_lb->nstage)
948 pme_lb->bBalance = FALSE;
950 if (DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
952 /* Unlock the DLB=auto, DLB is allowed to activate */
953 dd_dlb_unlock(cr->dd);
954 md_print_warn(cr, fp_log, "NOTE: DLB can now turn on, when beneficial\n");
956 /* We don't deactivate the tuning yet, since we will balance again
957 * after DLB gets turned on, if it does within PMETune_period.
959 continue_pme_loadbal(pme_lb, TRUE);
960 pme_lb->bTriggerOnDLB = TRUE;
961 pme_lb->step_rel_stop = step_rel + PMETunePeriod*ir->nstlist;
965 /* We're completely done with PME tuning */
966 pme_lb->bActive = FALSE;
969 if (DOMAINDECOMP(cr))
971 /* Set the cut-off limit to the final selected cut-off,
972 * so we don't have artificial DLB limits.
973 * This also ensures that we won't disable the currently
974 * optimal setting during a second round of PME balancing.
976 set_dd_dlb_max_cutoff(cr, fr->ic->rlist);
980 if (pme_lb->bBalance)
982 /* We might not have collected nstlist steps in cycles yet,
983 * since init_step might not be a multiple of nstlist,
984 * but the first data collected is skipped anyhow.
986 pme_load_balance(pme_lb, cr,
988 ir, state, pme_lb->cycles_c - cycles_prev,
989 fr->ic, fr->nbv, &fr->pmedata,
992 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
993 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
994 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
995 fr->rlist = fr->ic->rlist;
996 fr->rcoulomb = fr->ic->rcoulomb;
997 fr->rvdw = fr->ic->rvdw;
999 if (ir->eDispCorr != edispcNO)
1001 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1005 if (!pme_lb->bBalance &&
1006 (!pme_lb->bSepPMERanks || step_rel > pme_lb->step_rel_stop))
1008 /* We have just deactivated the balancing and we're not measuring PP/PME
1009 * imbalance during the first steps of the run: deactivate the tuning.
1011 pme_lb->bActive = FALSE;
1014 if (!(pme_lb->bActive) && DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
1016 /* Make sure DLB is allowed when we deactivate PME tuning */
1017 dd_dlb_unlock(cr->dd);
1018 md_print_warn(cr, fp_log, "NOTE: DLB can now turn on, when beneficial\n");
1021 *bPrinting = pme_lb->bBalance;
1024 /*! \brief Return product of the number of PME grid points in each dimension */
1025 static int pme_grid_points(const pme_setup_t *setup)
1027 return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
1030 /*! \brief Print one load-balancing setting */
1031 static void print_pme_loadbal_setting(FILE *fplog,
1033 const pme_setup_t *setup)
1036 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
1038 setup->rcut_coulomb, setup->rlist,
1039 setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
1040 setup->spacing, 1/setup->ewaldcoeff_q);
1043 /*! \brief Print all load-balancing settings */
1044 static void print_pme_loadbal_settings(pme_load_balancing_t *pme_lb,
1047 gmx_bool bNonBondedOnGPU)
1049 double pp_ratio, grid_ratio;
1050 real pp_ratio_temporary;
1052 pp_ratio_temporary = pme_lb->setup[pme_lb->cur].rlist / pme_lb->setup[0].rlist;
1053 pp_ratio = gmx::power3(pp_ratio_temporary);
1054 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
1055 (double)pme_grid_points(&pme_lb->setup[0]);
1057 fprintf(fplog, "\n");
1058 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
1059 fprintf(fplog, "\n");
1060 /* Here we only warn when the optimal setting is the last one */
1061 if (pme_lb->elimited != epmelblimNO &&
1062 pme_lb->cur == pme_loadbal_end(pme_lb)-1)
1064 fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
1065 pmelblim_str[pme_lb->elimited]);
1066 fprintf(fplog, " you might not have reached a good load balance.\n");
1067 if (pme_lb->elimited == epmelblimDD)
1069 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
1071 fprintf(fplog, "\n");
1073 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
1074 fprintf(fplog, " particle-particle PME\n");
1075 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
1076 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
1077 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
1078 fprintf(fplog, " cost-ratio %4.2f %4.2f\n",
1079 pp_ratio, grid_ratio);
1080 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
1082 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
1084 md_print_warn(cr, fplog,
1085 "NOTE: PME load balancing increased the non-bonded workload by more than 50%%.\n"
1086 " For better performance, use (more) PME ranks (mdrun -npme),\n"
1087 " or if you are beyond the scaling limit, use fewer total ranks (or nodes).\n");
1091 fprintf(fplog, "\n");
1095 void pme_loadbal_done(pme_load_balancing_t *pme_lb,
1098 gmx_bool bNonBondedOnGPU)
1100 if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
1102 print_pme_loadbal_settings(pme_lb, cr, fplog, bNonBondedOnGPU);
1105 /* TODO: Here we should free all pointers in pme_lb,
1106 * but as it contains pme data structures,
1107 * we need to first make pme.c free all data.