2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * \brief This file contains function definitions necessary for
38 * managing automatic load balance of PME calculations (Coulomb and
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_ewald
46 #include "pme-load-balancing.h"
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 */
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
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 */
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 int n
; /**< the count of setup as well as the allocation size */
153 pme_setup_t
*setup
; /**< the PME+cutoff setups */
154 int cur
; /**< the inex (in setup) of the current setup */
155 int fastest
; /**< index of the fastest setup up till now */
156 int lower_limit
; /**< don't go below this setup index */
157 int start
; /**< start of setup index range to consider in stage>0 */
158 int end
; /**< end of setup index range to consider in stage>0 */
159 int elimited
; /**< was the balancing limited, uses enum above */
160 int cutoff_scheme
; /**< Verlet or group cut-offs */
162 int stage
; /**< the current stage */
164 int cycles_n
; /**< step cycle counter cummulative count */
165 double cycles_c
; /**< step cycle counter cummulative cycles */
168 /* TODO The code in this file should call this getter, rather than
169 * read bActive anywhere */
170 bool pme_loadbal_is_active(const pme_load_balancing_t
*pme_lb
)
172 return pme_lb
!= nullptr && pme_lb
->bActive
;
175 void pme_loadbal_init(pme_load_balancing_t
**pme_lb_p
,
177 const gmx::MDLogger
&mdlog
,
178 const t_inputrec
&ir
,
180 const interaction_const_t
&ic
,
181 const NbnxnListParameters
&listParams
,
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
;
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");
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 */
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
);
222 snew(pme_lb
->setup
, pme_lb
->n
);
224 pme_lb
->rcut_vdw
= ic
.rvdw
;
225 pme_lb
->rcut_coulomb_start
= ir
.rcoulomb
;
228 pme_lb
->setup
[0].rcut_coulomb
= ic
.rcoulomb
;
229 pme_lb
->setup
[0].rlistOuter
= listParams
.rlistOuter
;
230 pme_lb
->setup
[0].rlistInner
= listParams
.rlistInner
;
231 pme_lb
->setup
[0].grid
[XX
] = ir
.nkx
;
232 pme_lb
->setup
[0].grid
[YY
] = ir
.nky
;
233 pme_lb
->setup
[0].grid
[ZZ
] = ir
.nkz
;
234 pme_lb
->setup
[0].ewaldcoeff_q
= ic
.ewaldcoeff_q
;
235 pme_lb
->setup
[0].ewaldcoeff_lj
= ic
.ewaldcoeff_lj
;
237 if (!pme_lb
->bSepPMERanks
)
239 GMX_RELEASE_ASSERT(pmedata
, "On ranks doing both PP and PME we need a valid pmedata object");
240 pme_lb
->setup
[0].pmedata
= pmedata
;
244 for (d
= 0; d
< DIM
; d
++)
246 sp
= norm(pme_lb
->box_start
[d
])/pme_lb
->setup
[0].grid
[d
];
252 pme_lb
->setup
[0].spacing
= spm
;
254 if (ir
.fourier_spacing
> 0)
256 pme_lb
->cut_spacing
= ir
.rcoulomb
/ir
.fourier_spacing
;
260 pme_lb
->cut_spacing
= ir
.rcoulomb
/pme_lb
->setup
[0].spacing
;
266 pme_lb
->lower_limit
= 0;
269 pme_lb
->elimited
= epmelblimNO
;
271 pme_lb
->cycles_n
= 0;
272 pme_lb
->cycles_c
= 0;
274 if (!wallcycle_have_counter())
276 GMX_LOG(mdlog
.warning
).asParagraph().appendText("NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use PME-PP balancing.");
279 /* Tune with GPUs and/or separate PME ranks.
280 * When running only on a CPU without PME ranks, PME tuning will only help
281 * with small numbers of atoms in the cut-off sphere.
283 pme_lb
->bActive
= (wallcycle_have_counter() && (bUseGPU
||
284 pme_lb
->bSepPMERanks
));
286 /* With GPUs and no separate PME ranks we can't measure the PP/PME
287 * imbalance, so we start balancing right away.
288 * Otherwise we only start balancing after we observe imbalance.
290 pme_lb
->bBalance
= (pme_lb
->bActive
&& (bUseGPU
&& !pme_lb
->bSepPMERanks
));
292 pme_lb
->step_rel_stop
= PMETunePeriod
*ir
.nstlist
;
294 /* Delay DD load balancing when GPUs are used */
295 if (pme_lb
->bActive
&& DOMAINDECOMP(cr
) && cr
->dd
->nnodes
> 1 && bUseGPU
)
297 /* Lock DLB=auto to off (does nothing when DLB=yes/no.
298 * With GPUs and separate PME nodes, we want to first
299 * do PME tuning without DLB, since DLB might limit
300 * the cut-off, which never improves performance.
301 * We allow for DLB + PME tuning after a first round of tuning.
304 if (dd_dlb_is_locked(cr
->dd
))
306 GMX_LOG(mdlog
.warning
).asParagraph().appendText("NOTE: DLB will not turn on during the first phase of PME tuning");
312 *bPrinting
= pme_lb
->bBalance
;
315 /*! \brief Try to increase the cutoff during load balancing */
316 static gmx_bool
pme_loadbal_increase_cutoff(pme_load_balancing_t
*pme_lb
,
318 const gmx_domdec_t
*dd
)
322 real tmpr_coulomb
, tmpr_vdw
;
326 /* Try to add a new setup with next larger cut-off to the list */
328 srenew(pme_lb
->setup
, pme_lb
->n
);
329 set
= &pme_lb
->setup
[pme_lb
->n
-1];
330 set
->pmedata
= nullptr;
332 NumPmeDomains numPmeDomains
= getNumPmeDomains(dd
);
337 /* Avoid infinite while loop, which can occur at the minimum grid size.
338 * Note that in practice load balancing will stop before this point.
339 * The factor 2.1 allows for the extreme case in which only grids
340 * of powers of 2 are allowed (the current code supports more grids).
350 clear_ivec(set
->grid
);
351 sp
= calcFftGrid(nullptr, pme_lb
->box_start
,
352 fac
*pme_lb
->setup
[pme_lb
->cur
].spacing
,
353 minimalPmeGridSize(pme_order
),
358 /* As here we can't easily check if one of the PME ranks
359 * uses threading, we do a conservative grid check.
360 * This means we can't use pme_order or less grid lines
361 * per PME rank along x, which is not a strong restriction.
363 grid_ok
= gmx_pme_check_restrictions(pme_order
,
364 set
->grid
[XX
], set
->grid
[YY
], set
->grid
[ZZ
],
369 while (sp
<= 1.001*pme_lb
->setup
[pme_lb
->cur
].spacing
|| !grid_ok
);
371 set
->rcut_coulomb
= pme_lb
->cut_spacing
*sp
;
372 if (set
->rcut_coulomb
< pme_lb
->rcut_coulomb_start
)
374 /* This is unlikely, but can happen when e.g. continuing from
375 * a checkpoint after equilibration where the box shrank a lot.
376 * We want to avoid rcoulomb getting smaller than rvdw
377 * and there might be more issues with decreasing rcoulomb.
379 set
->rcut_coulomb
= pme_lb
->rcut_coulomb_start
;
382 if (pme_lb
->cutoff_scheme
== ecutsVERLET
)
384 /* Never decrease the Coulomb and VdW list buffers */
385 set
->rlistOuter
= std::max(set
->rcut_coulomb
+ pme_lb
->rbufOuter_coulomb
,
386 pme_lb
->rcut_vdw
+ pme_lb
->rbufOuter_vdw
);
387 set
->rlistInner
= std::max(set
->rcut_coulomb
+ pme_lb
->rbufInner_coulomb
,
388 pme_lb
->rcut_vdw
+ pme_lb
->rbufInner_vdw
);
392 /* TODO Remove these lines and pme_lb->cutoff_scheme */
393 tmpr_coulomb
= set
->rcut_coulomb
+ pme_lb
->rbufOuter_coulomb
;
394 tmpr_vdw
= pme_lb
->rcut_vdw
+ pme_lb
->rbufOuter_vdw
;
395 /* Two (known) bugs with cutoff-scheme=group here:
396 * - This modification of rlist results in incorrect DD comunication.
397 * - We should set fr->bTwinRange = (fr->rlistlong > fr->rlist).
399 set
->rlistOuter
= std::min(tmpr_coulomb
, tmpr_vdw
);
400 set
->rlistInner
= set
->rlistOuter
;
404 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
405 set
->grid_efficiency
= 1;
406 for (d
= 0; d
< DIM
; d
++)
408 set
->grid_efficiency
*= (set
->grid
[d
]*sp
)/norm(pme_lb
->box_start
[d
]);
410 /* The Ewald coefficient is inversly proportional to the cut-off */
412 pme_lb
->setup
[0].ewaldcoeff_q
*pme_lb
->setup
[0].rcut_coulomb
/set
->rcut_coulomb
;
413 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
415 pme_lb
->setup
[0].ewaldcoeff_lj
*pme_lb
->setup
[0].rcut_coulomb
/set
->rcut_coulomb
;
422 fprintf(debug
, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
423 set
->grid
[XX
], set
->grid
[YY
], set
->grid
[ZZ
], set
->rcut_coulomb
);
428 /*! \brief Print the PME grid */
429 static void print_grid(FILE *fp_err
, FILE *fp_log
,
432 const pme_setup_t
*set
,
435 auto buf
= gmx::formatString("%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f",
437 set
->grid
[XX
], set
->grid
[YY
], set
->grid
[ZZ
], set
->rcut_coulomb
);
440 buf
+= gmx::formatString(": %.1f M-cycles", cycles
*1e-6);
442 if (fp_err
!= nullptr)
444 fprintf(fp_err
, "\r%s\n", buf
.c_str());
447 if (fp_log
!= nullptr)
449 fprintf(fp_log
, "%s\n", buf
.c_str());
453 /*! \brief Return the index of the last setup used in PME load balancing */
454 static int pme_loadbal_end(pme_load_balancing_t
*pme_lb
)
456 /* In the initial stage only n is set; end is not set yet */
467 /*! \brief Print descriptive string about what limits PME load balancing */
468 static void print_loadbal_limited(FILE *fp_err
, FILE *fp_log
,
470 pme_load_balancing_t
*pme_lb
)
472 auto buf
= gmx::formatString("step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
473 gmx::int64ToString(step
).c_str(),
474 pmelblim_str
[pme_lb
->elimited
],
475 pme_lb
->setup
[pme_loadbal_end(pme_lb
)-1].rcut_coulomb
);
476 if (fp_err
!= nullptr)
478 fprintf(fp_err
, "\r%s\n", buf
.c_str());
481 if (fp_log
!= nullptr)
483 fprintf(fp_log
, "%s\n", buf
.c_str());
487 /*! \brief Switch load balancing to stage 1
489 * In this stage, only reasonably fast setups are run again. */
490 static void switch_to_stage1(pme_load_balancing_t
*pme_lb
)
492 /* Increase start until we find a setup that is not slower than
493 * maxRelativeSlowdownAccepted times the fastest setup.
495 pme_lb
->start
= pme_lb
->lower_limit
;
496 while (pme_lb
->start
+ 1 < pme_lb
->n
&&
497 (pme_lb
->setup
[pme_lb
->start
].count
== 0 ||
498 pme_lb
->setup
[pme_lb
->start
].cycles
>
499 pme_lb
->setup
[pme_lb
->fastest
].cycles
*maxRelativeSlowdownAccepted
))
503 /* While increasing start, we might have skipped setups that we did not
504 * time during stage 0. We want to extend the range for stage 1 to include
505 * any skipped setups that lie between setups that were measured to be
506 * acceptably fast and too slow.
508 while (pme_lb
->start
> pme_lb
->lower_limit
&&
509 pme_lb
->setup
[pme_lb
->start
- 1].count
== 0)
514 /* Decrease end only with setups that we timed and that are slow. */
515 pme_lb
->end
= pme_lb
->n
;
516 if (pme_lb
->setup
[pme_lb
->end
- 1].count
> 0 &&
517 pme_lb
->setup
[pme_lb
->end
- 1].cycles
>
518 pme_lb
->setup
[pme_lb
->fastest
].cycles
*maxRelativeSlowdownAccepted
)
525 /* Next we want to choose setup pme_lb->end-1, but as we will decrease
526 * pme_lb->cur by one right after returning, we set cur to end.
528 pme_lb
->cur
= pme_lb
->end
;
531 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
533 * The adjustment is done to generate a different non-bonded PP and PME load.
534 * With separate PME ranks (PP and PME on different processes) or with
535 * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
536 * and changing the load will affect the load balance and performance.
537 * The total time for a set of integration steps is monitored and a range
538 * of grid/cut-off setups is scanned. After calling pme_load_balance many
539 * times and acquiring enough statistics, the best performing setup is chosen.
540 * Here we try to take into account fluctuations and changes due to external
541 * factors as well as DD load balancing.
544 pme_load_balance(pme_load_balancing_t
*pme_lb
,
548 const gmx::MDLogger
&mdlog
,
549 const t_inputrec
&ir
,
550 const t_state
&state
,
552 interaction_const_t
*ic
,
553 struct nonbonded_verlet_t
*nbv
,
554 struct gmx_pme_t
** pmedata
,
560 char buf
[STRLEN
], sbuf
[22];
565 gmx_sumd(1, &cycles
, cr
);
566 cycles
/= cr
->nnodes
;
569 set
= &pme_lb
->setup
[pme_lb
->cur
];
572 rtab
= ir
.rlist
+ ir
.tabext
;
574 if (set
->count
% 2 == 1)
576 /* Skip the first cycle, because the first step after a switch
577 * is much slower due to allocation and/or caching effects.
582 sprintf(buf
, "step %4s: ", gmx_step_str(step
, sbuf
));
583 print_grid(fp_err
, fp_log
, buf
, "timed with", set
, cycles
);
587 set
->cycles
= cycles
;
591 if (cycles
*maxFluctuationAccepted
< set
->cycles
&&
592 pme_lb
->stage
== pme_lb
->nstage
- 1)
594 /* The performance went up a lot (due to e.g. DD load balancing).
595 * Add a stage, keep the minima, but rescan all setups.
601 fprintf(debug
, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
602 "Increased the number stages to %d"
603 " and ignoring the previous performance\n",
604 set
->grid
[XX
], set
->grid
[YY
], set
->grid
[ZZ
],
605 set
->cycles
*1e-6, cycles
*1e-6, maxFluctuationAccepted
,
609 set
->cycles
= std::min(set
->cycles
, cycles
);
612 if (set
->cycles
< pme_lb
->setup
[pme_lb
->fastest
].cycles
)
614 pme_lb
->fastest
= pme_lb
->cur
;
616 if (DOMAINDECOMP(cr
))
618 /* We found a new fastest setting, ensure that with subsequent
619 * shorter cut-off's the dynamic load balancing does not make
620 * the use of the current cut-off impossible. This solution is
621 * a trade-off, as the PME load balancing and DD domain size
622 * load balancing can interact in complex ways.
623 * With the Verlet kernels, DD load imbalance will usually be
624 * mainly due to bonded interaction imbalance, which will often
625 * quickly push the domain boundaries beyond the limit for the
626 * optimal, PME load balanced, cut-off. But it could be that
627 * better overal performance can be obtained with a slightly
628 * shorter cut-off and better DD load balancing.
630 set_dd_dlb_max_cutoff(cr
, pme_lb
->setup
[pme_lb
->fastest
].rlistOuter
);
633 cycles_fast
= pme_lb
->setup
[pme_lb
->fastest
].cycles
;
635 /* Check in stage 0 if we should stop scanning grids.
636 * Stop when the time is more than maxRelativeSlowDownAccepted longer than the fastest.
638 if (pme_lb
->stage
== 0 && pme_lb
->cur
> 0 &&
639 cycles
> pme_lb
->setup
[pme_lb
->fastest
].cycles
*maxRelativeSlowdownAccepted
)
641 pme_lb
->n
= pme_lb
->cur
+ 1;
642 /* Done with scanning, go to stage 1 */
643 switch_to_stage1(pme_lb
);
646 if (pme_lb
->stage
== 0)
650 gridsize_start
= set
->grid
[XX
]*set
->grid
[YY
]*set
->grid
[ZZ
];
654 if (pme_lb
->cur
+1 < pme_lb
->n
)
656 /* We had already generated the next setup */
661 /* Find the next setup */
662 OK
= pme_loadbal_increase_cutoff(pme_lb
, ir
.pme_order
, cr
->dd
);
666 pme_lb
->elimited
= epmelblimPMEGRID
;
671 pme_lb
->setup
[pme_lb
->cur
+1].spacing
> c_maxSpacingScaling
*pme_lb
->setup
[0].spacing
)
674 pme_lb
->elimited
= epmelblimMAXSCALING
;
677 if (OK
&& ir
.ePBC
!= epbcNONE
)
679 OK
= (gmx::square(pme_lb
->setup
[pme_lb
->cur
+1].rlistOuter
)
680 <= max_cutoff2(ir
.ePBC
, state
.box
));
683 pme_lb
->elimited
= epmelblimBOX
;
691 if (DOMAINDECOMP(cr
))
693 OK
= change_dd_cutoff(cr
, state
,
694 pme_lb
->setup
[pme_lb
->cur
].rlistOuter
);
697 /* Failed: do not use this setup */
699 pme_lb
->elimited
= epmelblimDD
;
705 /* We hit the upper limit for the cut-off,
706 * the setup should not go further than cur.
708 pme_lb
->n
= pme_lb
->cur
+ 1;
709 print_loadbal_limited(fp_err
, fp_log
, step
, pme_lb
);
710 /* Switch to the next stage */
711 switch_to_stage1(pme_lb
);
715 !(pme_lb
->setup
[pme_lb
->cur
].grid
[XX
]*
716 pme_lb
->setup
[pme_lb
->cur
].grid
[YY
]*
717 pme_lb
->setup
[pme_lb
->cur
].grid
[ZZ
] <
718 gridsize_start
*gridpointsScaleFactor
720 pme_lb
->setup
[pme_lb
->cur
].grid_efficiency
<
721 pme_lb
->setup
[pme_lb
->cur
-1].grid_efficiency
*relativeEfficiencyFactor
));
724 if (pme_lb
->stage
> 0 && pme_lb
->end
== 1)
726 pme_lb
->cur
= pme_lb
->lower_limit
;
727 pme_lb
->stage
= pme_lb
->nstage
;
729 else if (pme_lb
->stage
> 0 && pme_lb
->end
> 1)
731 /* If stage = nstage-1:
732 * scan over all setups, rerunning only those setups
733 * which are not much slower than the fastest
736 * Note that we loop backward to minimize the risk of the cut-off
737 * getting limited by DD DLB, since the DLB cut-off limit is set
738 * to the fastest PME setup.
742 if (pme_lb
->cur
> pme_lb
->start
)
750 pme_lb
->cur
= pme_lb
->end
- 1;
753 while (pme_lb
->stage
== pme_lb
->nstage
- 1 &&
754 pme_lb
->setup
[pme_lb
->cur
].count
> 0 &&
755 pme_lb
->setup
[pme_lb
->cur
].cycles
> cycles_fast
*maxRelativeSlowdownAccepted
);
757 if (pme_lb
->stage
== pme_lb
->nstage
)
759 /* We are done optimizing, use the fastest setup we found */
760 pme_lb
->cur
= pme_lb
->fastest
;
764 if (DOMAINDECOMP(cr
) && pme_lb
->stage
> 0)
766 OK
= change_dd_cutoff(cr
, state
, pme_lb
->setup
[pme_lb
->cur
].rlistOuter
);
769 /* For some reason the chosen cut-off is incompatible with DD.
770 * We should continue scanning a more limited range of cut-off's.
772 if (pme_lb
->cur
> 1 && pme_lb
->stage
== pme_lb
->nstage
)
774 /* stage=nstage says we're finished, but we should continue
775 * balancing, so we set back stage which was just incremented.
779 if (pme_lb
->cur
<= pme_lb
->fastest
)
781 /* This should not happen, as we set limits on the DLB bounds.
782 * But we implement a complete failsafe solution anyhow.
784 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
785 "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
);
786 pme_lb
->fastest
= pme_lb
->lower_limit
;
787 pme_lb
->start
= pme_lb
->lower_limit
;
789 /* Limit the range to below the current cut-off, scan from start */
790 pme_lb
->end
= pme_lb
->cur
;
791 pme_lb
->cur
= pme_lb
->start
;
792 pme_lb
->elimited
= epmelblimDD
;
793 print_loadbal_limited(fp_err
, fp_log
, step
, pme_lb
);
797 /* Change the Coulomb cut-off and the PME grid */
799 set
= &pme_lb
->setup
[pme_lb
->cur
];
801 NbnxnListParameters
*listParams
= nbv
->listParams
.get();
803 ic
->rcoulomb
= set
->rcut_coulomb
;
804 listParams
->rlistOuter
= set
->rlistOuter
;
805 listParams
->rlistInner
= set
->rlistInner
;
806 ic
->ewaldcoeff_q
= set
->ewaldcoeff_q
;
807 /* TODO: centralize the code that sets the potentials shifts */
808 if (ic
->coulomb_modifier
== eintmodPOTSHIFT
)
810 GMX_RELEASE_ASSERT(ic
->rcoulomb
!= 0, "Cutoff radius cannot be zero");
811 ic
->sh_ewald
= std::erfc(ic
->ewaldcoeff_q
*ic
->rcoulomb
) / ic
->rcoulomb
;
813 if (EVDW_PME(ic
->vdwtype
))
815 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
816 ic
->rvdw
= set
->rcut_coulomb
;
817 ic
->ewaldcoeff_lj
= set
->ewaldcoeff_lj
;
818 if (ic
->vdw_modifier
== eintmodPOTSHIFT
)
822 ic
->dispersion_shift
.cpot
= -1.0/gmx::power6(static_cast<double>(ic
->rvdw
));
823 ic
->repulsion_shift
.cpot
= -1.0/gmx::power12(static_cast<double>(ic
->rvdw
));
824 ic
->sh_invrc6
= -ic
->dispersion_shift
.cpot
;
825 crc2
= gmx::square(ic
->ewaldcoeff_lj
*ic
->rvdw
);
826 ic
->sh_lj_ewald
= (std::exp(-crc2
)*(1 + crc2
+ 0.5*crc2
*crc2
) - 1)/gmx::power6(ic
->rvdw
);
830 /* We always re-initialize the tables whether they are used or not */
831 init_interaction_const_tables(nullptr, ic
, rtab
);
833 nbnxn_gpu_pme_loadbal_update_param(nbv
, ic
, listParams
);
835 if (!pme_lb
->bSepPMERanks
)
838 * CPU PME keeps a list of allocated pmedata's, that's why pme_lb->setup[pme_lb->cur].pmedata is not always nullptr.
839 * GPU PME, however, currently needs the gmx_pme_reinit always called on load balancing
840 * (pme_gpu_reinit might be not sufficiently decoupled from gmx_pme_init).
841 * This can lead to a lot of reallocations for PME GPU.
842 * Would be nicer if the allocated grid list was hidden within a single pmedata structure.
844 if ((pme_lb
->setup
[pme_lb
->cur
].pmedata
== nullptr) || pme_gpu_task_enabled(pme_lb
->setup
[pme_lb
->cur
].pmedata
))
846 /* Generate a new PME data structure,
847 * copying part of the old pointers.
849 gmx_pme_reinit(&set
->pmedata
,
850 cr
, pme_lb
->setup
[0].pmedata
, &ir
,
851 set
->grid
, set
->ewaldcoeff_q
, set
->ewaldcoeff_lj
);
853 *pmedata
= set
->pmedata
;
857 /* Tell our PME-only rank to switch grid */
858 gmx_pme_send_switchgrid(cr
, set
->grid
, set
->ewaldcoeff_q
, set
->ewaldcoeff_lj
);
863 print_grid(nullptr, debug
, "", "switched to", set
, -1);
866 if (pme_lb
->stage
== pme_lb
->nstage
)
868 print_grid(fp_err
, fp_log
, "", "optimal", set
, -1);
872 /*! \brief Prepare for another round of PME load balancing
874 * \param[in,out] pme_lb Pointer to PME load balancing struct
875 * \param[in] bDlbUnlocked TRUE is DLB was locked and is now unlocked
877 * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
878 * the PP/PME balance might change and re-balancing can improve performance.
879 * This function adds 2 stages and adjusts the considered setup range.
881 static void continue_pme_loadbal(pme_load_balancing_t
*pme_lb
,
882 gmx_bool bDlbUnlocked
)
884 /* Add 2 tuning stages, keep the detected end of the setup range */
886 if (bDlbUnlocked
&& pme_lb
->bSepPMERanks
)
888 /* With separate PME ranks, DLB should always lower the PP load and
889 * can only increase the PME load (more communication and imbalance),
890 * so we only need to scan longer cut-off's.
892 pme_lb
->lower_limit
= pme_lb
->cur
;
894 pme_lb
->start
= pme_lb
->lower_limit
;
897 void pme_loadbal_do(pme_load_balancing_t
*pme_lb
,
901 const gmx::MDLogger
&mdlog
,
902 const t_inputrec
&ir
,
904 const t_state
&state
,
905 gmx_wallcycle_t wcycle
,
913 assert(pme_lb
!= nullptr);
915 if (!pme_lb
->bActive
)
920 n_prev
= pme_lb
->cycles_n
;
921 cycles_prev
= pme_lb
->cycles_c
;
922 wallcycle_get(wcycle
, ewcSTEP
, &pme_lb
->cycles_n
, &pme_lb
->cycles_c
);
924 /* Before the first step we haven't done any steps yet.
925 * Also handle cases where ir.init_step % ir.nstlist != 0.
927 if (pme_lb
->cycles_n
< ir
.nstlist
)
931 /* Sanity check, we expect nstlist cycle counts */
932 if (pme_lb
->cycles_n
- n_prev
!= ir
.nstlist
)
934 /* We could return here, but it's safer to issue an error and quit */
935 gmx_incons("pme_loadbal_do called at an interval != nstlist");
938 /* PME grid + cut-off optimization with GPUs or PME ranks */
939 if (!pme_lb
->bBalance
&& pme_lb
->bSepPMERanks
)
941 if (pme_lb
->bTriggerOnDLB
)
943 pme_lb
->bBalance
= dd_dlb_is_on(cr
->dd
);
945 /* We should ignore the first timing to avoid timing allocation
946 * overhead. And since the PME load balancing is called just
947 * before DD repartitioning, the ratio returned by dd_pme_f_ratio
948 * is not over the last nstlist steps, but the nstlist steps before
949 * that. So the first useful ratio is available at step_rel=3*nstlist.
951 else if (step_rel
>= 3*ir
.nstlist
)
953 if (DDMASTER(cr
->dd
))
955 /* If PME rank load is too high, start tuning */
957 (dd_pme_f_ratio(cr
->dd
) >= loadBalanceTriggerFactor
);
959 dd_bcast(cr
->dd
, sizeof(gmx_bool
), &pme_lb
->bBalance
);
962 pme_lb
->bActive
= (pme_lb
->bBalance
||
963 step_rel
<= pme_lb
->step_rel_stop
);
966 /* The location in the code of this balancing termination is strange.
967 * You would expect to have it after the call to pme_load_balance()
968 * below, since there pme_lb->stage is updated.
969 * But when terminating directly after deciding on and selecting the
970 * optimal setup, DLB will turn on right away if it was locked before.
971 * This might be due to PME reinitialization. So we check stage here
972 * to allow for another nstlist steps with DLB locked to stabilize
975 if (pme_lb
->bBalance
&& pme_lb
->stage
== pme_lb
->nstage
)
977 pme_lb
->bBalance
= FALSE
;
979 if (DOMAINDECOMP(cr
) && dd_dlb_is_locked(cr
->dd
))
981 /* Unlock the DLB=auto, DLB is allowed to activate */
982 dd_dlb_unlock(cr
->dd
);
983 GMX_LOG(mdlog
.warning
).asParagraph().appendText("NOTE: DLB can now turn on, when beneficial");
985 /* We don't deactivate the tuning yet, since we will balance again
986 * after DLB gets turned on, if it does within PMETune_period.
988 continue_pme_loadbal(pme_lb
, TRUE
);
989 pme_lb
->bTriggerOnDLB
= TRUE
;
990 pme_lb
->step_rel_stop
= step_rel
+ PMETunePeriod
*ir
.nstlist
;
994 /* We're completely done with PME tuning */
995 pme_lb
->bActive
= FALSE
;
998 if (DOMAINDECOMP(cr
))
1000 /* Set the cut-off limit to the final selected cut-off,
1001 * so we don't have artificial DLB limits.
1002 * This also ensures that we won't disable the currently
1003 * optimal setting during a second round of PME balancing.
1005 set_dd_dlb_max_cutoff(cr
, fr
->nbv
->listParams
->rlistOuter
);
1009 if (pme_lb
->bBalance
)
1011 /* We might not have collected nstlist steps in cycles yet,
1012 * since init_step might not be a multiple of nstlist,
1013 * but the first data collected is skipped anyhow.
1015 pme_load_balance(pme_lb
, cr
,
1016 fp_err
, fp_log
, mdlog
,
1017 ir
, state
, pme_lb
->cycles_c
- cycles_prev
,
1018 fr
->ic
, fr
->nbv
, &fr
->pmedata
,
1021 /* Update deprecated rlist in forcerec to stay in sync with fr->nbv */
1022 fr
->rlist
= fr
->nbv
->listParams
->rlistOuter
;
1024 if (ir
.eDispCorr
!= edispcNO
)
1026 calc_enervirdiff(nullptr, ir
.eDispCorr
, fr
);
1030 if (!pme_lb
->bBalance
&&
1031 (!pme_lb
->bSepPMERanks
|| step_rel
> pme_lb
->step_rel_stop
))
1033 /* We have just deactivated the balancing and we're not measuring PP/PME
1034 * imbalance during the first steps of the run: deactivate the tuning.
1036 pme_lb
->bActive
= FALSE
;
1039 if (!(pme_lb
->bActive
) && DOMAINDECOMP(cr
) && dd_dlb_is_locked(cr
->dd
))
1041 /* Make sure DLB is allowed when we deactivate PME tuning */
1042 dd_dlb_unlock(cr
->dd
);
1043 GMX_LOG(mdlog
.warning
).asParagraph().appendText("NOTE: DLB can now turn on, when beneficial");
1046 *bPrinting
= pme_lb
->bBalance
;
1049 /*! \brief Return product of the number of PME grid points in each dimension */
1050 static int pme_grid_points(const pme_setup_t
*setup
)
1052 return setup
->grid
[XX
]*setup
->grid
[YY
]*setup
->grid
[ZZ
];
1055 /*! \brief Print one load-balancing setting */
1056 static void print_pme_loadbal_setting(FILE *fplog
,
1058 const pme_setup_t
*setup
)
1061 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
1063 setup
->rcut_coulomb
, setup
->rlistInner
,
1064 setup
->grid
[XX
], setup
->grid
[YY
], setup
->grid
[ZZ
],
1065 setup
->spacing
, 1/setup
->ewaldcoeff_q
);
1068 /*! \brief Print all load-balancing settings */
1069 static void print_pme_loadbal_settings(pme_load_balancing_t
*pme_lb
,
1071 const gmx::MDLogger
&mdlog
,
1072 gmx_bool bNonBondedOnGPU
)
1074 double pp_ratio
, grid_ratio
;
1075 real pp_ratio_temporary
;
1077 pp_ratio_temporary
= pme_lb
->setup
[pme_lb
->cur
].rlistInner
/ pme_lb
->setup
[0].rlistInner
;
1078 pp_ratio
= gmx::power3(pp_ratio_temporary
);
1079 grid_ratio
= pme_grid_points(&pme_lb
->setup
[pme_lb
->cur
])/
1080 static_cast<double>(pme_grid_points(&pme_lb
->setup
[0]));
1082 fprintf(fplog
, "\n");
1083 fprintf(fplog
, " P P - P M E L O A D B A L A N C I N G\n");
1084 fprintf(fplog
, "\n");
1085 /* Here we only warn when the optimal setting is the last one */
1086 if (pme_lb
->elimited
!= epmelblimNO
&&
1087 pme_lb
->cur
== pme_loadbal_end(pme_lb
)-1)
1089 fprintf(fplog
, " NOTE: The PP/PME load balancing was limited by the %s,\n",
1090 pmelblim_str
[pme_lb
->elimited
]);
1091 fprintf(fplog
, " you might not have reached a good load balance.\n");
1092 if (pme_lb
->elimited
== epmelblimDD
)
1094 fprintf(fplog
, " Try different mdrun -dd settings or lower the -dds value.\n");
1096 fprintf(fplog
, "\n");
1098 fprintf(fplog
, " PP/PME load balancing changed the cut-off and PME settings:\n");
1099 fprintf(fplog
, " particle-particle PME\n");
1100 fprintf(fplog
, " rcoulomb rlist grid spacing 1/beta\n");
1101 print_pme_loadbal_setting(fplog
, "initial", &pme_lb
->setup
[0]);
1102 print_pme_loadbal_setting(fplog
, "final", &pme_lb
->setup
[pme_lb
->cur
]);
1103 fprintf(fplog
, " cost-ratio %4.2f %4.2f\n",
1104 pp_ratio
, grid_ratio
);
1105 fprintf(fplog
, " (note that these numbers concern only part of the total PP and PME load)\n");
1107 if (pp_ratio
> 1.5 && !bNonBondedOnGPU
)
1109 GMX_LOG(mdlog
.warning
).asParagraph().appendText(
1110 "NOTE: PME load balancing increased the non-bonded workload by more than 50%.\n"
1111 " For better performance, use (more) PME ranks (mdrun -npme),\n"
1112 " or if you are beyond the scaling limit, use fewer total ranks (or nodes).");
1116 fprintf(fplog
, "\n");
1120 void pme_loadbal_done(pme_load_balancing_t
*pme_lb
,
1122 const gmx::MDLogger
&mdlog
,
1123 gmx_bool bNonBondedOnGPU
)
1125 if (fplog
!= nullptr && (pme_lb
->cur
> 0 || pme_lb
->elimited
!= epmelblimNO
))
1127 print_pme_loadbal_settings(pme_lb
, fplog
, mdlog
, bNonBondedOnGPU
);
1130 /* TODO: Here we should free all pointers in pme_lb,
1131 * but as it contains pme data structures,
1132 * we need to first make pme.c free all data.