Simple union-find structure
[alexxy/gromacs.git] / src / gromacs / ewald / pme-load-balancing.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief This file contains function definitions necessary for
38  * managing automatic load balance of PME calculations (Coulomb and
39  * LJ).
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_ewald
43  */
44 #include "gmxpre.h"
45
46 #include "pme-load-balancing.h"
47
48 #include "config.h"
49
50 #include <assert.h>
51
52 #include <cmath>
53
54 #include <algorithm>
55
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.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
81 #include "pme-internal.h"
82
83 /*! \brief Parameters and settings for one PP-PME setup */
84 struct pme_setup_t {
85     real              rcut_coulomb;    /**< Coulomb cut-off                              */
86     real              rlistOuter;      /**< cut-off for the outer pair-list              */
87     real              rlistInner;      /**< cut-off for the inner pair-list              */
88     real              spacing;         /**< (largest) PME grid spacing                   */
89     ivec              grid;            /**< the PME grid dimensions                      */
90     real              grid_efficiency; /**< ineffiency factor for non-uniform grids <= 1 */
91     real              ewaldcoeff_q;    /**< Electrostatic Ewald coefficient            */
92     real              ewaldcoeff_lj;   /**< LJ Ewald coefficient, only for the call to send_switchgrid */
93     struct gmx_pme_t *pmedata;         /**< the data structure used in the PME code      */
94     int               count;           /**< number of times this setup has been timed    */
95     double            cycles;          /**< the fastest time for this setup in cycles    */
96 };
97
98 /*! \brief After 50 nstlist periods of not observing imbalance: never tune PME */
99 const int  PMETunePeriod = 50;
100 /*! \brief Trigger PME load balancing at more than 5% PME overload */
101 const real loadBalanceTriggerFactor = 1.05;
102 /*! \brief In the initial scan, step by grids that are at least a factor 0.8 coarser */
103 const real gridScaleFactor = 0.8;
104 /*! \brief In the initial scan, try to skip grids with uneven x/y/z spacing,
105  * checking if the "efficiency" is more than 5% worse than the previous grid.
106  */
107 const real relativeEfficiencyFactor = 1.05;
108 /*! \brief Rerun until a run is 12% slower setups than the fastest run so far */
109 const real maxRelativeSlowdownAccepted = 1.12;
110 /*! \brief If setups get more than 2% faster, do another round to avoid
111  * choosing a slower setup due to acceleration or fluctuations.
112  */
113 const real maxFluctuationAccepted = 1.02;
114
115 /*! \brief Enumeration whose values describe the effect limiting the load balancing */
116 enum epmelb {
117     epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimPMEGRID, epmelblimNR
118 };
119
120 /*! \brief Descriptive strings matching ::epmelb */
121 const char *pmelblim_str[epmelblimNR] =
122 { "no", "box size", "domain decompostion", "PME grid restriction" };
123
124 struct pme_load_balancing_t {
125     gmx_bool     bSepPMERanks;       /**< do we have separate PME ranks? */
126     gmx_bool     bActive;            /**< is PME tuning active? */
127     gmx_int64_t  step_rel_stop;      /**< stop the tuning after this value of step_rel */
128     gmx_bool     bTriggerOnDLB;      /**< trigger balancing only on DD DLB */
129     gmx_bool     bBalance;           /**< are we in the balancing phase, i.e. trying different setups? */
130     int          nstage;             /**< the current maximum number of stages */
131
132     real         cut_spacing;        /**< the minimum cutoff / PME grid spacing ratio */
133     real         rcut_vdw;           /**< Vdw cutoff (does not change) */
134     real         rcut_coulomb_start; /**< Initial electrostatics cutoff */
135     real         rbufOuter_coulomb;  /**< the outer pairlist buffer size */
136     real         rbufOuter_vdw;      /**< the outer pairlist buffer size */
137     real         rbufInner_coulomb;  /**< the inner pairlist buffer size */
138     real         rbufInner_vdw;      /**< the inner pairlist buffer size */
139     matrix       box_start;          /**< the initial simulation box */
140     int          n;                  /**< the count of setup as well as the allocation size */
141     pme_setup_t *setup;              /**< the PME+cutoff setups */
142     int          cur;                /**< the inex (in setup) of the current setup */
143     int          fastest;            /**< index of the fastest setup up till now */
144     int          lower_limit;        /**< don't go below this setup index */
145     int          start;              /**< start of setup index range to consider in stage>0 */
146     int          end;                /**< end   of setup index range to consider in stage>0 */
147     int          elimited;           /**< was the balancing limited, uses enum above */
148     int          cutoff_scheme;      /**< Verlet or group cut-offs */
149
150     int          stage;              /**< the current stage */
151
152     int          cycles_n;           /**< step cycle counter cummulative count */
153     double       cycles_c;           /**< step cycle counter cummulative cycles */
154 };
155
156 /* TODO The code in this file should call this getter, rather than
157  * read bActive anywhere */
158 bool pme_loadbal_is_active(const pme_load_balancing_t *pme_lb)
159 {
160     return pme_lb != nullptr && pme_lb->bActive;
161 }
162
163 void pme_loadbal_init(pme_load_balancing_t     **pme_lb_p,
164                       t_commrec                 *cr,
165                       const gmx::MDLogger       &mdlog,
166                       const t_inputrec          *ir,
167                       matrix                     box,
168                       const interaction_const_t *ic,
169                       const NbnxnListParameters *listParams,
170                       gmx_pme_t                 *pmedata,
171                       gmx_bool                   bUseGPU,
172                       gmx_bool                  *bPrinting)
173 {
174     pme_load_balancing_t *pme_lb;
175     real                  spm, sp;
176     int                   d;
177
178     // Note that we don't (yet) support PME load balancing with LJ-PME only.
179     GMX_RELEASE_ASSERT(EEL_PME(ir->coulombtype), "pme_loadbal_init called without PME electrostatics");
180     // To avoid complexity, we require a single cut-off with PME for q+LJ.
181     // This is checked by grompp, but it doesn't hurt to check again.
182     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");
183
184     snew(pme_lb, 1);
185
186     pme_lb->bSepPMERanks      = !(cr->duty & DUTY_PME);
187
188     /* Initially we turn on balancing directly on based on PP/PME imbalance */
189     pme_lb->bTriggerOnDLB     = FALSE;
190
191     /* Any number of stages >= 2 is supported */
192     pme_lb->nstage            = 2;
193
194     pme_lb->cutoff_scheme     = ir->cutoff_scheme;
195
196     pme_lb->rbufOuter_coulomb = listParams->rlistOuter - ic->rcoulomb;
197     pme_lb->rbufOuter_vdw     = listParams->rlistOuter - ic->rvdw;
198     pme_lb->rbufInner_coulomb = listParams->rlistInner - ic->rcoulomb;
199     pme_lb->rbufInner_vdw     = listParams->rlistInner - ic->rvdw;
200
201     copy_mat(box, pme_lb->box_start);
202     if (ir->ePBC == epbcXY && ir->nwall == 2)
203     {
204         svmul(ir->wall_ewald_zfac, pme_lb->box_start[ZZ], pme_lb->box_start[ZZ]);
205     }
206
207     pme_lb->n = 1;
208     snew(pme_lb->setup, pme_lb->n);
209
210     pme_lb->rcut_vdw                 = ic->rvdw;
211     pme_lb->rcut_coulomb_start       = ir->rcoulomb;
212
213     pme_lb->cur                      = 0;
214     pme_lb->setup[0].rcut_coulomb    = ic->rcoulomb;
215     pme_lb->setup[0].rlistOuter      = listParams->rlistOuter;
216     pme_lb->setup[0].rlistInner      = listParams->rlistInner;
217     pme_lb->setup[0].grid[XX]        = ir->nkx;
218     pme_lb->setup[0].grid[YY]        = ir->nky;
219     pme_lb->setup[0].grid[ZZ]        = ir->nkz;
220     pme_lb->setup[0].ewaldcoeff_q    = ic->ewaldcoeff_q;
221     pme_lb->setup[0].ewaldcoeff_lj   = ic->ewaldcoeff_lj;
222
223     pme_lb->setup[0].pmedata         = pmedata;
224
225     spm = 0;
226     for (d = 0; d < DIM; d++)
227     {
228         sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
229         if (sp > spm)
230         {
231             spm = sp;
232         }
233     }
234     pme_lb->setup[0].spacing = spm;
235
236     if (ir->fourier_spacing > 0)
237     {
238         pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
239     }
240     else
241     {
242         pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
243     }
244
245     pme_lb->stage = 0;
246
247     pme_lb->fastest     = 0;
248     pme_lb->lower_limit = 0;
249     pme_lb->start       = 0;
250     pme_lb->end         = 0;
251     pme_lb->elimited    = epmelblimNO;
252
253     pme_lb->cycles_n = 0;
254     pme_lb->cycles_c = 0;
255
256     if (!wallcycle_have_counter())
257     {
258         GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: Cycle counters unsupported or not enabled in kernel. Cannot use PME-PP balancing.");
259     }
260
261     /* Tune with GPUs and/or separate PME ranks.
262      * When running only on a CPU without PME ranks, PME tuning will only help
263      * with small numbers of atoms in the cut-off sphere.
264      */
265     pme_lb->bActive  = (wallcycle_have_counter() && (bUseGPU ||
266                                                      pme_lb->bSepPMERanks));
267
268     /* With GPUs and no separate PME ranks we can't measure the PP/PME
269      * imbalance, so we start balancing right away.
270      * Otherwise we only start balancing after we observe imbalance.
271      */
272     pme_lb->bBalance = (pme_lb->bActive && (bUseGPU && !pme_lb->bSepPMERanks));
273
274     pme_lb->step_rel_stop = PMETunePeriod*ir->nstlist;
275
276     /* Delay DD load balancing when GPUs are used */
277     if (pme_lb->bActive && DOMAINDECOMP(cr) && cr->dd->nnodes > 1 && bUseGPU)
278     {
279         /* Lock DLB=auto to off (does nothing when DLB=yes/no.
280          * With GPUs and separate PME nodes, we want to first
281          * do PME tuning without DLB, since DLB might limit
282          * the cut-off, which never improves performance.
283          * We allow for DLB + PME tuning after a first round of tuning.
284          */
285         dd_dlb_lock(cr->dd);
286         if (dd_dlb_is_locked(cr->dd))
287         {
288             GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: DLB will not turn on during the first phase of PME tuning");
289         }
290     }
291
292     *pme_lb_p = pme_lb;
293
294     *bPrinting = pme_lb->bBalance;
295 }
296
297 /*! \brief Try to increase the cutoff during load balancing */
298 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t *pme_lb,
299                                             int                   pme_order,
300                                             const gmx_domdec_t   *dd)
301 {
302     pme_setup_t *set;
303     int          npmeranks_x, npmeranks_y;
304     real         fac, sp;
305     real         tmpr_coulomb, tmpr_vdw;
306     int          d;
307     bool         grid_ok;
308
309     /* Try to add a new setup with next larger cut-off to the list */
310     pme_lb->n++;
311     srenew(pme_lb->setup, pme_lb->n);
312     set          = &pme_lb->setup[pme_lb->n-1];
313     set->pmedata = nullptr;
314
315     get_pme_nnodes(dd, &npmeranks_x, &npmeranks_y);
316
317     fac = 1;
318     do
319     {
320         /* Avoid infinite while loop, which can occur at the minimum grid size.
321          * Note that in practice load balancing will stop before this point.
322          * The factor 2.1 allows for the extreme case in which only grids
323          * of powers of 2 are allowed (the current code supports more grids).
324          */
325         if (fac > 2.1)
326         {
327             pme_lb->n--;
328
329             return FALSE;
330         }
331
332         fac *= 1.01;
333         clear_ivec(set->grid);
334         sp = calcFftGrid(nullptr, pme_lb->box_start,
335                          fac*pme_lb->setup[pme_lb->cur].spacing,
336                          minimalPmeGridSize(pme_order),
337                          &set->grid[XX],
338                          &set->grid[YY],
339                          &set->grid[ZZ]);
340
341         /* As here we can't easily check if one of the PME ranks
342          * uses threading, we do a conservative grid check.
343          * This means we can't use pme_order or less grid lines
344          * per PME rank along x, which is not a strong restriction.
345          */
346         grid_ok = gmx_pme_check_restrictions(pme_order,
347                                              set->grid[XX], set->grid[YY], set->grid[ZZ],
348                                              npmeranks_x,
349                                              true,
350                                              false);
351     }
352     while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
353
354     set->rcut_coulomb = pme_lb->cut_spacing*sp;
355     if (set->rcut_coulomb < pme_lb->rcut_coulomb_start)
356     {
357         /* This is unlikely, but can happen when e.g. continuing from
358          * a checkpoint after equilibration where the box shrank a lot.
359          * We want to avoid rcoulomb getting smaller than rvdw
360          * and there might be more issues with decreasing rcoulomb.
361          */
362         set->rcut_coulomb = pme_lb->rcut_coulomb_start;
363     }
364
365     if (pme_lb->cutoff_scheme == ecutsVERLET)
366     {
367         /* Never decrease the Coulomb and VdW list buffers */
368         set->rlistOuter  = std::max(set->rcut_coulomb + pme_lb->rbufOuter_coulomb,
369                                     pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw);
370         set->rlistInner  = std::max(set->rcut_coulomb + pme_lb->rbufInner_coulomb,
371                                     pme_lb->rcut_vdw + pme_lb->rbufInner_vdw);
372     }
373     else
374     {
375         tmpr_coulomb     = set->rcut_coulomb + pme_lb->rbufOuter_coulomb;
376         tmpr_vdw         = pme_lb->rcut_vdw + pme_lb->rbufOuter_vdw;
377         set->rlistOuter  = std::min(tmpr_coulomb, tmpr_vdw);
378         set->rlistInner  = set->rlistOuter;
379     }
380
381     set->spacing         = sp;
382     /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
383     set->grid_efficiency = 1;
384     for (d = 0; d < DIM; d++)
385     {
386         set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
387     }
388     /* The Ewald coefficient is inversly proportional to the cut-off */
389     set->ewaldcoeff_q =
390         pme_lb->setup[0].ewaldcoeff_q*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
391     /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
392     set->ewaldcoeff_lj =
393         pme_lb->setup[0].ewaldcoeff_lj*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
394
395     set->count   = 0;
396     set->cycles  = 0;
397
398     if (debug)
399     {
400         fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
401                 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
402     }
403     return TRUE;
404 }
405
406 /*! \brief Print the PME grid */
407 static void print_grid(FILE *fp_err, FILE *fp_log,
408                        const char *pre,
409                        const char *desc,
410                        const pme_setup_t *set,
411                        double cycles)
412 {
413     char buf[STRLEN], buft[STRLEN];
414
415     if (cycles >= 0)
416     {
417         sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
418     }
419     else
420     {
421         buft[0] = '\0';
422     }
423     sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
424             pre,
425             desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
426             buft);
427     if (fp_err != nullptr)
428     {
429         fprintf(fp_err, "\r%s\n", buf);
430         fflush(fp_err);
431     }
432     if (fp_log != nullptr)
433     {
434         fprintf(fp_log, "%s\n", buf);
435     }
436 }
437
438 /*! \brief Return the index of the last setup used in PME load balancing */
439 static int pme_loadbal_end(pme_load_balancing_t *pme_lb)
440 {
441     /* In the initial stage only n is set; end is not set yet */
442     if (pme_lb->end > 0)
443     {
444         return pme_lb->end;
445     }
446     else
447     {
448         return pme_lb->n;
449     }
450 }
451
452 /*! \brief Print descriptive string about what limits PME load balancing */
453 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
454                                   gmx_int64_t step,
455                                   pme_load_balancing_t *pme_lb)
456 {
457     char buf[STRLEN], sbuf[22];
458
459     sprintf(buf, "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
460             gmx_step_str(step, sbuf),
461             pmelblim_str[pme_lb->elimited],
462             pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
463     if (fp_err != nullptr)
464     {
465         fprintf(fp_err, "\r%s\n", buf);
466         fflush(fp_err);
467     }
468     if (fp_log != nullptr)
469     {
470         fprintf(fp_log, "%s\n", buf);
471     }
472 }
473
474 /*! \brief Switch load balancing to stage 1
475  *
476  * In this stage, only reasonably fast setups are run again. */
477 static void switch_to_stage1(pme_load_balancing_t *pme_lb)
478 {
479     /* Increase start until we find a setup that is not slower than
480      * maxRelativeSlowdownAccepted times the fastest setup.
481      */
482     pme_lb->start = pme_lb->lower_limit;
483     while (pme_lb->start + 1 < pme_lb->n &&
484            (pme_lb->setup[pme_lb->start].count == 0 ||
485             pme_lb->setup[pme_lb->start].cycles >
486             pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted))
487     {
488         pme_lb->start++;
489     }
490     /* While increasing start, we might have skipped setups that we did not
491      * time during stage 0. We want to extend the range for stage 1 to include
492      * any skipped setups that lie between setups that were measured to be
493      * acceptably fast and too slow.
494      */
495     while (pme_lb->start > pme_lb->lower_limit &&
496            pme_lb->setup[pme_lb->start - 1].count == 0)
497     {
498         pme_lb->start--;
499     }
500
501     /* Decrease end only with setups that we timed and that are slow. */
502     pme_lb->end = pme_lb->n;
503     if (pme_lb->setup[pme_lb->end - 1].count > 0 &&
504         pme_lb->setup[pme_lb->end - 1].cycles >
505         pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted)
506     {
507         pme_lb->end--;
508     }
509
510     pme_lb->stage = 1;
511
512     /* Next we want to choose setup pme_lb->end-1, but as we will decrease
513      * pme_lb->cur by one right after returning, we set cur to end.
514      */
515     pme_lb->cur = pme_lb->end;
516 }
517
518 /*! \brief Process the timings and try to adjust the PME grid and Coulomb cut-off
519  *
520  * The adjustment is done to generate a different non-bonded PP and PME load.
521  * With separate PME ranks (PP and PME on different processes) or with
522  * a GPU (PP on GPU, PME on CPU), PP and PME run on different resources
523  * and changing the load will affect the load balance and performance.
524  * The total time for a set of integration steps is monitored and a range
525  * of grid/cut-off setups is scanned. After calling pme_load_balance many
526  * times and acquiring enough statistics, the best performing setup is chosen.
527  * Here we try to take into account fluctuations and changes due to external
528  * factors as well as DD load balancing.
529  */
530 static void
531 pme_load_balance(pme_load_balancing_t      *pme_lb,
532                  t_commrec                 *cr,
533                  FILE                      *fp_err,
534                  FILE                      *fp_log,
535                  const gmx::MDLogger       &mdlog,
536                  const t_inputrec          *ir,
537                  t_state                   *state,
538                  double                     cycles,
539                  interaction_const_t       *ic,
540                  struct nonbonded_verlet_t *nbv,
541                  struct gmx_pme_t **        pmedata,
542                  gmx_int64_t                step)
543 {
544     gmx_bool     OK;
545     pme_setup_t *set;
546     double       cycles_fast;
547     char         buf[STRLEN], sbuf[22];
548     real         rtab;
549
550     if (PAR(cr))
551     {
552         gmx_sumd(1, &cycles, cr);
553         cycles /= cr->nnodes;
554     }
555
556     set = &pme_lb->setup[pme_lb->cur];
557     set->count++;
558
559     rtab = ir->rlist + ir->tabext;
560
561     if (set->count % 2 == 1)
562     {
563         /* Skip the first cycle, because the first step after a switch
564          * is much slower due to allocation and/or caching effects.
565          */
566         return;
567     }
568
569     sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
570     print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
571
572     if (set->count <= 2)
573     {
574         set->cycles = cycles;
575     }
576     else
577     {
578         if (cycles*maxFluctuationAccepted < set->cycles &&
579             pme_lb->stage == pme_lb->nstage - 1)
580         {
581             /* The performance went up a lot (due to e.g. DD load balancing).
582              * Add a stage, keep the minima, but rescan all setups.
583              */
584             pme_lb->nstage++;
585
586             if (debug)
587             {
588                 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
589                         "Increased the number stages to %d"
590                         " and ignoring the previous performance\n",
591                         set->grid[XX], set->grid[YY], set->grid[ZZ],
592                         set->cycles*1e-6, cycles*1e-6, maxFluctuationAccepted,
593                         pme_lb->nstage);
594             }
595         }
596         set->cycles = std::min(set->cycles, cycles);
597     }
598
599     if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
600     {
601         pme_lb->fastest = pme_lb->cur;
602
603         if (DOMAINDECOMP(cr))
604         {
605             /* We found a new fastest setting, ensure that with subsequent
606              * shorter cut-off's the dynamic load balancing does not make
607              * the use of the current cut-off impossible. This solution is
608              * a trade-off, as the PME load balancing and DD domain size
609              * load balancing can interact in complex ways.
610              * With the Verlet kernels, DD load imbalance will usually be
611              * mainly due to bonded interaction imbalance, which will often
612              * quickly push the domain boundaries beyond the limit for the
613              * optimal, PME load balanced, cut-off. But it could be that
614              * better overal performance can be obtained with a slightly
615              * shorter cut-off and better DD load balancing.
616              */
617             set_dd_dlb_max_cutoff(cr, pme_lb->setup[pme_lb->fastest].rlistOuter);
618         }
619     }
620     cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
621
622     /* Check in stage 0 if we should stop scanning grids.
623      * Stop when the time is more than maxRelativeSlowDownAccepted longer than the fastest.
624      */
625     if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
626         cycles > pme_lb->setup[pme_lb->fastest].cycles*maxRelativeSlowdownAccepted)
627     {
628         pme_lb->n = pme_lb->cur + 1;
629         /* Done with scanning, go to stage 1 */
630         switch_to_stage1(pme_lb);
631     }
632
633     if (pme_lb->stage == 0)
634     {
635         int gridsize_start;
636
637         gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
638
639         do
640         {
641             if (pme_lb->cur+1 < pme_lb->n)
642             {
643                 /* We had already generated the next setup */
644                 OK = TRUE;
645             }
646             else
647             {
648                 /* Find the next setup */
649                 OK = pme_loadbal_increase_cutoff(pme_lb, ir->pme_order, cr->dd);
650
651                 if (!OK)
652                 {
653                     pme_lb->elimited = epmelblimPMEGRID;
654                 }
655             }
656
657             if (OK && ir->ePBC != epbcNONE)
658             {
659                 OK = (gmx::square(pme_lb->setup[pme_lb->cur+1].rlistOuter)
660                       <= max_cutoff2(ir->ePBC, state->box));
661                 if (!OK)
662                 {
663                     pme_lb->elimited = epmelblimBOX;
664                 }
665             }
666
667             if (OK)
668             {
669                 pme_lb->cur++;
670
671                 if (DOMAINDECOMP(cr))
672                 {
673                     OK = change_dd_cutoff(cr, state, ir,
674                                           pme_lb->setup[pme_lb->cur].rlistOuter);
675                     if (!OK)
676                     {
677                         /* Failed: do not use this setup */
678                         pme_lb->cur--;
679                         pme_lb->elimited = epmelblimDD;
680                     }
681                 }
682             }
683             if (!OK)
684             {
685                 /* We hit the upper limit for the cut-off,
686                  * the setup should not go further than cur.
687                  */
688                 pme_lb->n = pme_lb->cur + 1;
689                 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
690                 /* Switch to the next stage */
691                 switch_to_stage1(pme_lb);
692             }
693         }
694         while (OK &&
695                !(pme_lb->setup[pme_lb->cur].grid[XX]*
696                  pme_lb->setup[pme_lb->cur].grid[YY]*
697                  pme_lb->setup[pme_lb->cur].grid[ZZ] <
698                  gridsize_start*gridScaleFactor
699                  &&
700                  pme_lb->setup[pme_lb->cur].grid_efficiency <
701                  pme_lb->setup[pme_lb->cur-1].grid_efficiency*relativeEfficiencyFactor));
702     }
703
704     if (pme_lb->stage > 0 && pme_lb->end == 1)
705     {
706         pme_lb->cur   = pme_lb->lower_limit;
707         pme_lb->stage = pme_lb->nstage;
708     }
709     else if (pme_lb->stage > 0 && pme_lb->end > 1)
710     {
711         /* If stage = nstage-1:
712          *   scan over all setups, rerunning only those setups
713          *   which are not much slower than the fastest
714          * else:
715          *   use the next setup
716          * Note that we loop backward to minimize the risk of the cut-off
717          * getting limited by DD DLB, since the DLB cut-off limit is set
718          * to the fastest PME setup.
719          */
720         do
721         {
722             if (pme_lb->cur > pme_lb->start)
723             {
724                 pme_lb->cur--;
725             }
726             else
727             {
728                 pme_lb->stage++;
729
730                 pme_lb->cur = pme_lb->end - 1;
731             }
732         }
733         while (pme_lb->stage == pme_lb->nstage - 1 &&
734                pme_lb->setup[pme_lb->cur].count > 0 &&
735                pme_lb->setup[pme_lb->cur].cycles > cycles_fast*maxRelativeSlowdownAccepted);
736
737         if (pme_lb->stage == pme_lb->nstage)
738         {
739             /* We are done optimizing, use the fastest setup we found */
740             pme_lb->cur = pme_lb->fastest;
741         }
742     }
743
744     if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
745     {
746         OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlistOuter);
747         if (!OK)
748         {
749             /* For some reason the chosen cut-off is incompatible with DD.
750              * We should continue scanning a more limited range of cut-off's.
751              */
752             if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
753             {
754                 /* stage=nstage says we're finished, but we should continue
755                  * balancing, so we set back stage which was just incremented.
756                  */
757                 pme_lb->stage--;
758             }
759             if (pme_lb->cur <= pme_lb->fastest)
760             {
761                 /* This should not happen, as we set limits on the DLB bounds.
762                  * But we implement a complete failsafe solution anyhow.
763                  */
764                 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
765                         "The fastest PP/PME load balancing setting (cutoff %.3f nm) is no longer available due to DD DLB or box size limitations", pme_lb->fastest);
766                 pme_lb->fastest = pme_lb->lower_limit;
767                 pme_lb->start   = pme_lb->lower_limit;
768             }
769             /* Limit the range to below the current cut-off, scan from start */
770             pme_lb->end         = pme_lb->cur;
771             pme_lb->cur         = pme_lb->start;
772             pme_lb->elimited    = epmelblimDD;
773             print_loadbal_limited(fp_err, fp_log, step, pme_lb);
774         }
775     }
776
777     /* Change the Coulomb cut-off and the PME grid */
778
779     set = &pme_lb->setup[pme_lb->cur];
780
781     NbnxnListParameters *listParams = nbv->listParams.get();
782
783     ic->rcoulomb           = set->rcut_coulomb;
784     listParams->rlistOuter = set->rlistOuter;
785     listParams->rlistInner = set->rlistInner;
786     ic->ewaldcoeff_q       = set->ewaldcoeff_q;
787     /* TODO: centralize the code that sets the potentials shifts */
788     if (ic->coulomb_modifier == eintmodPOTSHIFT)
789     {
790         ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb);
791     }
792     if (EVDW_PME(ic->vdwtype))
793     {
794         /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
795         ic->rvdw            = set->rcut_coulomb;
796         ic->ewaldcoeff_lj   = set->ewaldcoeff_lj;
797         if (ic->vdw_modifier == eintmodPOTSHIFT)
798         {
799             real       crc2;
800
801             ic->dispersion_shift.cpot = -1.0/gmx::power6(static_cast<double>(ic->rvdw));
802             ic->repulsion_shift.cpot  = -1.0/gmx::power12(static_cast<double>(ic->rvdw));
803             ic->sh_invrc6             = -ic->dispersion_shift.cpot;
804             crc2                      = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
805             ic->sh_lj_ewald           = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
806         }
807     }
808
809     /* We always re-initialize the tables whether they are used or not */
810     init_interaction_const_tables(nullptr, ic, rtab);
811
812     nbnxn_gpu_pme_loadbal_update_param(nbv, ic, listParams);
813
814     /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
815      * also sharing texture references. To keep the code simple, we don't
816      * treat texture references as shared resources, but this means that
817      * the coulomb_tab texture ref will get updated by multiple threads.
818      * Hence, to ensure that the non-bonded kernels don't start before all
819      * texture binding operations are finished, we need to wait for all ranks
820      * to arrive here before continuing.
821      *
822      * Note that we could omit this barrier if GPUs are not shared (or
823      * texture objects are used), but as this is initialization code, there
824      * is not point in complicating things.
825      */
826 #if GMX_THREAD_MPI
827     if (PAR(cr) && use_GPU(nbv))
828     {
829         gmx_barrier(cr);
830     }
831 #endif  /* GMX_THREAD_MPI */
832
833     if (!pme_lb->bSepPMERanks)
834     {
835         if (pme_lb->setup[pme_lb->cur].pmedata == nullptr)
836         {
837             /* Generate a new PME data structure,
838              * copying part of the old pointers.
839              */
840             gmx_pme_reinit(&set->pmedata,
841                            cr, pme_lb->setup[0].pmedata, ir,
842                            set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
843         }
844         *pmedata = set->pmedata;
845     }
846     else
847     {
848         /* Tell our PME-only rank to switch grid */
849         gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
850     }
851
852     if (debug)
853     {
854         print_grid(nullptr, debug, "", "switched to", set, -1);
855     }
856
857     if (pme_lb->stage == pme_lb->nstage)
858     {
859         print_grid(fp_err, fp_log, "", "optimal", set, -1);
860     }
861 }
862
863 /*! \brief Prepare for another round of PME load balancing
864  *
865  * \param[in,out] pme_lb       Pointer to PME load balancing struct
866  * \param[in]     bDlbUnlocked TRUE is DLB was locked and is now unlocked
867  *
868  * If the conditions (e.g. DLB off/on, CPU/GPU throttling etc.) changed,
869  * the PP/PME balance might change and re-balancing can improve performance.
870  * This function adds 2 stages and adjusts the considered setup range.
871  */
872 static void continue_pme_loadbal(pme_load_balancing_t *pme_lb,
873                                  gmx_bool              bDlbUnlocked)
874 {
875     /* Add 2 tuning stages, keep the detected end of the setup range */
876     pme_lb->nstage          += 2;
877     if (bDlbUnlocked && pme_lb->bSepPMERanks)
878     {
879         /* With separate PME ranks, DLB should always lower the PP load and
880          * can only increase the PME load (more communication and imbalance),
881          * so we only need to scan longer cut-off's.
882          */
883         pme_lb->lower_limit  = pme_lb->cur;
884     }
885     pme_lb->start            = pme_lb->lower_limit;
886 }
887
888 void pme_loadbal_do(pme_load_balancing_t *pme_lb,
889                     t_commrec            *cr,
890                     FILE                 *fp_err,
891                     FILE                 *fp_log,
892                     const gmx::MDLogger  &mdlog,
893                     const t_inputrec     *ir,
894                     t_forcerec           *fr,
895                     t_state              *state,
896                     gmx_wallcycle_t       wcycle,
897                     gmx_int64_t           step,
898                     gmx_int64_t           step_rel,
899                     gmx_bool             *bPrinting)
900 {
901     int    n_prev;
902     double cycles_prev;
903
904     assert(pme_lb != NULL);
905
906     if (!pme_lb->bActive)
907     {
908         return;
909     }
910
911     n_prev      = pme_lb->cycles_n;
912     cycles_prev = pme_lb->cycles_c;
913     wallcycle_get(wcycle, ewcSTEP, &pme_lb->cycles_n, &pme_lb->cycles_c);
914
915     /* Before the first step we haven't done any steps yet.
916      * Also handle cases where ir->init_step % ir->nstlist != 0.
917      */
918     if (pme_lb->cycles_n < ir->nstlist)
919     {
920         return;
921     }
922     /* Sanity check, we expect nstlist cycle counts */
923     if (pme_lb->cycles_n - n_prev != ir->nstlist)
924     {
925         /* We could return here, but it's safer to issue an error and quit */
926         gmx_incons("pme_loadbal_do called at an interval != nstlist");
927     }
928
929     /* PME grid + cut-off optimization with GPUs or PME ranks */
930     if (!pme_lb->bBalance && pme_lb->bSepPMERanks)
931     {
932         if (pme_lb->bTriggerOnDLB)
933         {
934             pme_lb->bBalance = dd_dlb_is_on(cr->dd);
935         }
936         /* We should ignore the first timing to avoid timing allocation
937          * overhead. And since the PME load balancing is called just
938          * before DD repartitioning, the ratio returned by dd_pme_f_ratio
939          * is not over the last nstlist steps, but the nstlist steps before
940          * that. So the first useful ratio is available at step_rel=3*nstlist.
941          */
942         else if (step_rel >= 3*ir->nstlist)
943         {
944             if (DDMASTER(cr->dd))
945             {
946                 /* If PME rank load is too high, start tuning */
947                 pme_lb->bBalance =
948                     (dd_pme_f_ratio(cr->dd) >= loadBalanceTriggerFactor);
949             }
950             dd_bcast(cr->dd, sizeof(gmx_bool), &pme_lb->bBalance);
951         }
952
953         pme_lb->bActive = (pme_lb->bBalance ||
954                            step_rel <= pme_lb->step_rel_stop);
955     }
956
957     /* The location in the code of this balancing termination is strange.
958      * You would expect to have it after the call to pme_load_balance()
959      * below, since there pme_lb->stage is updated.
960      * But when terminating directly after deciding on and selecting the
961      * optimal setup, DLB will turn on right away if it was locked before.
962      * This might be due to PME reinitialization. So we check stage here
963      * to allow for another nstlist steps with DLB locked to stabilize
964      * the performance.
965      */
966     if (pme_lb->bBalance && pme_lb->stage == pme_lb->nstage)
967     {
968         pme_lb->bBalance = FALSE;
969
970         if (DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
971         {
972             /* Unlock the DLB=auto, DLB is allowed to activate */
973             dd_dlb_unlock(cr->dd);
974             GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: DLB can now turn on, when beneficial");
975
976             /* We don't deactivate the tuning yet, since we will balance again
977              * after DLB gets turned on, if it does within PMETune_period.
978              */
979             continue_pme_loadbal(pme_lb, TRUE);
980             pme_lb->bTriggerOnDLB = TRUE;
981             pme_lb->step_rel_stop = step_rel + PMETunePeriod*ir->nstlist;
982         }
983         else
984         {
985             /* We're completely done with PME tuning */
986             pme_lb->bActive = FALSE;
987         }
988
989         if (DOMAINDECOMP(cr))
990         {
991             /* Set the cut-off limit to the final selected cut-off,
992              * so we don't have artificial DLB limits.
993              * This also ensures that we won't disable the currently
994              * optimal setting during a second round of PME balancing.
995              */
996             set_dd_dlb_max_cutoff(cr, fr->nbv->listParams->rlistOuter);
997         }
998     }
999
1000     if (pme_lb->bBalance)
1001     {
1002         /* We might not have collected nstlist steps in cycles yet,
1003          * since init_step might not be a multiple of nstlist,
1004          * but the first data collected is skipped anyhow.
1005          */
1006         pme_load_balance(pme_lb, cr,
1007                          fp_err, fp_log, mdlog,
1008                          ir, state, pme_lb->cycles_c - cycles_prev,
1009                          fr->ic, fr->nbv, &fr->pmedata,
1010                          step);
1011
1012         /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1013         fr->ewaldcoeff_q  = fr->ic->ewaldcoeff_q;
1014         fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1015         fr->rlist         = fr->nbv->listParams->rlistOuter;
1016         fr->rcoulomb      = fr->ic->rcoulomb;
1017         fr->rvdw          = fr->ic->rvdw;
1018
1019         if (ir->eDispCorr != edispcNO)
1020         {
1021             calc_enervirdiff(nullptr, ir->eDispCorr, fr);
1022         }
1023     }
1024
1025     if (!pme_lb->bBalance &&
1026         (!pme_lb->bSepPMERanks || step_rel > pme_lb->step_rel_stop))
1027     {
1028         /* We have just deactivated the balancing and we're not measuring PP/PME
1029          * imbalance during the first steps of the run: deactivate the tuning.
1030          */
1031         pme_lb->bActive = FALSE;
1032     }
1033
1034     if (!(pme_lb->bActive) && DOMAINDECOMP(cr) && dd_dlb_is_locked(cr->dd))
1035     {
1036         /* Make sure DLB is allowed when we deactivate PME tuning */
1037         dd_dlb_unlock(cr->dd);
1038         GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: DLB can now turn on, when beneficial");
1039     }
1040
1041     *bPrinting = pme_lb->bBalance;
1042 }
1043
1044 /*! \brief Return product of the number of PME grid points in each dimension */
1045 static int pme_grid_points(const pme_setup_t *setup)
1046 {
1047     return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
1048 }
1049
1050 /*! \brief Print one load-balancing setting */
1051 static void print_pme_loadbal_setting(FILE              *fplog,
1052                                       const char        *name,
1053                                       const pme_setup_t *setup)
1054 {
1055     fprintf(fplog,
1056             "   %-7s %6.3f nm %6.3f nm     %3d %3d %3d   %5.3f nm  %5.3f nm\n",
1057             name,
1058             setup->rcut_coulomb, setup->rlistInner,
1059             setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
1060             setup->spacing, 1/setup->ewaldcoeff_q);
1061 }
1062
1063 /*! \brief Print all load-balancing settings */
1064 static void print_pme_loadbal_settings(pme_load_balancing_t *pme_lb,
1065                                        FILE                 *fplog,
1066                                        const gmx::MDLogger  &mdlog,
1067                                        gmx_bool              bNonBondedOnGPU)
1068 {
1069     double     pp_ratio, grid_ratio;
1070     real       pp_ratio_temporary;
1071
1072     pp_ratio_temporary = pme_lb->setup[pme_lb->cur].rlistInner / pme_lb->setup[0].rlistInner;
1073     pp_ratio           = gmx::power3(pp_ratio_temporary);
1074     grid_ratio         = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
1075         (double)pme_grid_points(&pme_lb->setup[0]);
1076
1077     fprintf(fplog, "\n");
1078     fprintf(fplog, "       P P   -   P M E   L O A D   B A L A N C I N G\n");
1079     fprintf(fplog, "\n");
1080     /* Here we only warn when the optimal setting is the last one */
1081     if (pme_lb->elimited != epmelblimNO &&
1082         pme_lb->cur == pme_loadbal_end(pme_lb)-1)
1083     {
1084         fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
1085                 pmelblim_str[pme_lb->elimited]);
1086         fprintf(fplog, "       you might not have reached a good load balance.\n");
1087         if (pme_lb->elimited == epmelblimDD)
1088         {
1089             fprintf(fplog, "       Try different mdrun -dd settings or lower the -dds value.\n");
1090         }
1091         fprintf(fplog, "\n");
1092     }
1093     fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
1094     fprintf(fplog, "           particle-particle                    PME\n");
1095     fprintf(fplog, "            rcoulomb  rlist            grid      spacing   1/beta\n");
1096     print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
1097     print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
1098     fprintf(fplog, " cost-ratio           %4.2f             %4.2f\n",
1099             pp_ratio, grid_ratio);
1100     fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
1101
1102     if (pp_ratio > 1.5 && !bNonBondedOnGPU)
1103     {
1104         GMX_LOG(mdlog.warning).asParagraph().appendText(
1105                 "NOTE: PME load balancing increased the non-bonded workload by more than 50%.\n"
1106                 "      For better performance, use (more) PME ranks (mdrun -npme),\n"
1107                 "      or if you are beyond the scaling limit, use fewer total ranks (or nodes).");
1108     }
1109     else
1110     {
1111         fprintf(fplog, "\n");
1112     }
1113 }
1114
1115 void pme_loadbal_done(pme_load_balancing_t *pme_lb,
1116                       FILE                 *fplog,
1117                       const gmx::MDLogger  &mdlog,
1118                       gmx_bool              bNonBondedOnGPU)
1119 {
1120     if (fplog != nullptr && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
1121     {
1122         print_pme_loadbal_settings(pme_lb, fplog, mdlog, bNonBondedOnGPU);
1123     }
1124
1125     /* TODO: Here we should free all pointers in pme_lb,
1126      * but as it contains pme data structures,
1127      * we need to first make pme.c free all data.
1128      */
1129 }