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