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