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