PME load balancing now works with LJ-PME
[alexxy/gromacs.git] / src / programs / mdrun / pme_loadbal.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, 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 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "smalloc.h"
40 #include "network.h"
41 #include "calcgrid.h"
42 #include "pme.h"
43 #include "vec.h"
44 #include "domdec.h"
45 #include "nbnxn_cuda_data_mgmt.h"
46 #include "force.h"
47 #include "macros.h"
48 #include "md_logging.h"
49 #include "pme_loadbal.h"
50
51 /* Parameters and setting for one PP-PME setup */
52 typedef struct {
53     real      rcut_coulomb;    /* Coulomb cut-off                              */
54     real      rlist;           /* pair-list cut-off                            */
55     real      rlistlong;       /* LR pair-list cut-off                         */
56     int       nstcalclr;       /* frequency of evaluating long-range forces for group scheme */
57     real      spacing;         /* (largest) PME grid spacing                   */
58     ivec      grid;            /* the PME grid dimensions                      */
59     real      grid_efficiency; /* ineffiency factor for non-uniform grids <= 1 */
60     real      ewaldcoeff_q;    /* Electrostatic Ewald coefficient            */
61     real      ewaldcoeff_lj;   /* LJ Ewald coefficient, only for the call to send_switchgrid */
62     gmx_pme_t pmedata;         /* the data structure used in the PME code      */
63     int       count;           /* number of times this setup has been timed    */
64     double    cycles;          /* the fastest time for this setup in cycles    */
65 } pme_setup_t;
66
67 /* In the initial scan, step by grids that are at least a factor 0.8 coarser */
68 #define PME_LB_GRID_SCALE_FAC  0.8
69 /* In the initial scan, try to skip grids with uneven x/y/z spacing,
70  * checking if the "efficiency" is more than 5% worse than the previous grid.
71  */
72 #define PME_LB_GRID_EFFICIENCY_REL_FAC  1.05
73 /* Rerun up till 12% slower setups than the fastest up till now */
74 #define PME_LB_SLOW_FAC  1.12
75 /* If setups get more than 2% faster, do another round to avoid
76  * choosing a slower setup due to acceleration or fluctuations.
77  */
78 #define PME_LB_ACCEL_TOL 1.02
79
80 enum {
81     epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimPMEGRID, epmelblimNR
82 };
83
84 const char *pmelblim_str[epmelblimNR] =
85 { "no", "box size", "domain decompostion", "PME grid restriction" };
86
87 struct pme_load_balancing {
88     int          nstage;             /* the current maximum number of stages */
89
90     real         cut_spacing;        /* the minimum cutoff / PME grid spacing ratio */
91     real         rcut_vdw;           /* Vdw cutoff (does not change) */
92     real         rcut_coulomb_start; /* Initial electrostatics cutoff */
93     int          nstcalclr_start;    /* Initial electrostatics cutoff */
94     real         rbuf_coulomb;       /* the pairlist buffer size */
95     real         rbuf_vdw;           /* the pairlist buffer size */
96     matrix       box_start;          /* the initial simulation box */
97     int          n;                  /* the count of setup as well as the allocation size */
98     pme_setup_t *setup;              /* the PME+cutoff setups */
99     int          cur;                /* the current setup */
100     int          fastest;            /* fastest setup up till now */
101     int          start;              /* start of setup range to consider in stage>0 */
102     int          end;                /* end   of setup range to consider in stage>0 */
103     int          elimited;           /* was the balancing limited, uses enum above */
104     int          cutoff_scheme;      /* Verlet or group cut-offs */
105
106     int          stage;              /* the current stage */
107 };
108
109 void pme_loadbal_init(pme_load_balancing_t *pme_lb_p,
110                       const t_inputrec *ir, matrix box,
111                       const interaction_const_t *ic,
112                       gmx_pme_t pmedata)
113 {
114     pme_load_balancing_t pme_lb;
115     real                 spm, sp;
116     int                  d;
117
118     snew(pme_lb, 1);
119
120     /* Any number of stages >= 2 is supported */
121     pme_lb->nstage   = 2;
122
123     pme_lb->cutoff_scheme = ir->cutoff_scheme;
124
125     if (pme_lb->cutoff_scheme == ecutsVERLET)
126     {
127         pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
128         pme_lb->rbuf_vdw     = pme_lb->rbuf_coulomb;
129     }
130     else
131     {
132         if (ic->rcoulomb > ic->rlist)
133         {
134             pme_lb->rbuf_coulomb = ic->rlistlong - ic->rcoulomb;
135         }
136         else
137         {
138             pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
139         }
140         if (ic->rvdw > ic->rlist)
141         {
142             pme_lb->rbuf_vdw = ic->rlistlong - ic->rvdw;
143         }
144         else
145         {
146             pme_lb->rbuf_vdw = ic->rlist - ic->rvdw;
147         }
148     }
149
150     copy_mat(box, pme_lb->box_start);
151     if (ir->ePBC == epbcXY && ir->nwall == 2)
152     {
153         svmul(ir->wall_ewald_zfac, pme_lb->box_start[ZZ], pme_lb->box_start[ZZ]);
154     }
155
156     pme_lb->n = 1;
157     snew(pme_lb->setup, pme_lb->n);
158
159     pme_lb->rcut_vdw                 = ic->rvdw;
160     pme_lb->rcut_coulomb_start       = ir->rcoulomb;
161     pme_lb->nstcalclr_start          = ir->nstcalclr;
162
163     pme_lb->cur                      = 0;
164     pme_lb->setup[0].rcut_coulomb    = ic->rcoulomb;
165     pme_lb->setup[0].rlist           = ic->rlist;
166     pme_lb->setup[0].rlistlong       = ic->rlistlong;
167     pme_lb->setup[0].nstcalclr       = ir->nstcalclr;
168     pme_lb->setup[0].grid[XX]        = ir->nkx;
169     pme_lb->setup[0].grid[YY]        = ir->nky;
170     pme_lb->setup[0].grid[ZZ]        = ir->nkz;
171     pme_lb->setup[0].ewaldcoeff_q    = ic->ewaldcoeff_q;
172     pme_lb->setup[0].ewaldcoeff_lj   = ic->ewaldcoeff_lj;
173
174     pme_lb->setup[0].pmedata  = pmedata;
175
176     spm = 0;
177     for (d = 0; d < DIM; d++)
178     {
179         sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
180         if (sp > spm)
181         {
182             spm = sp;
183         }
184     }
185     pme_lb->setup[0].spacing = spm;
186
187     if (ir->fourier_spacing > 0)
188     {
189         pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
190     }
191     else
192     {
193         pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
194     }
195
196     pme_lb->stage = 0;
197
198     pme_lb->fastest  = 0;
199     pme_lb->start    = 0;
200     pme_lb->end      = 0;
201     pme_lb->elimited = epmelblimNO;
202
203     *pme_lb_p = pme_lb;
204 }
205
206 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t  pme_lb,
207                                             int                   pme_order,
208                                             const gmx_domdec_t   *dd)
209 {
210     pme_setup_t *set;
211     int          npmenodes_x, npmenodes_y;
212     real         fac, sp;
213     real         tmpr_coulomb, tmpr_vdw;
214     int          d;
215     gmx_bool     grid_ok;
216
217     /* Try to add a new setup with next larger cut-off to the list */
218     pme_lb->n++;
219     srenew(pme_lb->setup, pme_lb->n);
220     set          = &pme_lb->setup[pme_lb->n-1];
221     set->pmedata = NULL;
222
223     get_pme_nnodes(dd, &npmenodes_x, &npmenodes_y);
224
225     fac = 1;
226     do
227     {
228         /* Avoid infinite while loop, which can occur at the minimum grid size.
229          * Note that in practice load balancing will stop before this point.
230          * The factor 2.1 allows for the extreme case in which only grids
231          * of powers of 2 are allowed (the current code supports more grids).
232          */
233         if (fac > 2.1)
234         {
235             pme_lb->n--;
236
237             return FALSE;
238         }
239
240         fac *= 1.01;
241         clear_ivec(set->grid);
242         sp = calc_grid(NULL, pme_lb->box_start,
243                        fac*pme_lb->setup[pme_lb->cur].spacing,
244                        &set->grid[XX],
245                        &set->grid[YY],
246                        &set->grid[ZZ]);
247
248         /* As here we can't easily check if one of the PME nodes
249          * uses threading, we do a conservative grid check.
250          * This means we can't use pme_order or less grid lines
251          * per PME node along x, which is not a strong restriction.
252          */
253         gmx_pme_check_restrictions(pme_order,
254                                    set->grid[XX], set->grid[YY], set->grid[ZZ],
255                                    npmenodes_x, npmenodes_y,
256                                    TRUE,
257                                    FALSE,
258                                    &grid_ok);
259     }
260     while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
261
262     set->rcut_coulomb = pme_lb->cut_spacing*sp;
263
264     if (pme_lb->cutoff_scheme == ecutsVERLET)
265     {
266         set->rlist        = set->rcut_coulomb + pme_lb->rbuf_coulomb;
267         /* We dont use LR lists with Verlet, but this avoids if-statements in further checks */
268         set->rlistlong    = set->rlist;
269     }
270     else
271     {
272         tmpr_coulomb          = set->rcut_coulomb + pme_lb->rbuf_coulomb;
273         tmpr_vdw              = pme_lb->rcut_vdw + pme_lb->rbuf_vdw;
274         set->rlist            = min(tmpr_coulomb, tmpr_vdw);
275         set->rlistlong        = max(tmpr_coulomb, tmpr_vdw);
276
277         /* Set the long-range update frequency */
278         if (set->rlist == set->rlistlong)
279         {
280             /* No long-range interactions if the short-/long-range cutoffs are identical */
281             set->nstcalclr = 0;
282         }
283         else if (pme_lb->nstcalclr_start == 0 || pme_lb->nstcalclr_start == 1)
284         {
285             /* We were not doing long-range before, but now we are since rlist!=rlistlong */
286             set->nstcalclr = 1;
287         }
288         else
289         {
290             /* We were already doing long-range interactions from the start */
291             if (pme_lb->rcut_vdw > pme_lb->rcut_coulomb_start)
292             {
293                 /* We were originally doing long-range VdW-only interactions.
294                  * If rvdw is still longer than rcoulomb we keep the original nstcalclr,
295                  * but if the coulomb cutoff has become longer we should update the long-range
296                  * part every step.
297                  */
298                 set->nstcalclr = (tmpr_vdw > tmpr_coulomb) ? pme_lb->nstcalclr_start : 1;
299             }
300             else
301             {
302                 /* We were not doing any long-range interaction from the start,
303                  * since it is not possible to do twin-range coulomb for the PME interaction.
304                  */
305                 set->nstcalclr = 1;
306             }
307         }
308     }
309
310     set->spacing      = sp;
311     /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
312     set->grid_efficiency = 1;
313     for (d = 0; d < DIM; d++)
314     {
315         set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
316     }
317     /* The Ewald coefficient is inversly proportional to the cut-off */
318     set->ewaldcoeff_q =
319         pme_lb->setup[0].ewaldcoeff_q*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
320     /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
321     set->ewaldcoeff_lj =
322         pme_lb->setup[0].ewaldcoeff_lj*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
323
324     set->count   = 0;
325     set->cycles  = 0;
326
327     if (debug)
328     {
329         fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
330                 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
331     }
332     return TRUE;
333 }
334
335 static void print_grid(FILE *fp_err, FILE *fp_log,
336                        const char *pre,
337                        const char *desc,
338                        const pme_setup_t *set,
339                        double cycles)
340 {
341     char buf[STRLEN], buft[STRLEN];
342
343     if (cycles >= 0)
344     {
345         sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
346     }
347     else
348     {
349         buft[0] = '\0';
350     }
351     sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
352             pre,
353             desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
354             buft);
355     if (fp_err != NULL)
356     {
357         fprintf(fp_err, "\r%s\n", buf);
358     }
359     if (fp_log != NULL)
360     {
361         fprintf(fp_log, "%s\n", buf);
362     }
363 }
364
365 static int pme_loadbal_end(pme_load_balancing_t pme_lb)
366 {
367     /* In the initial stage only n is set; end is not set yet */
368     if (pme_lb->end > 0)
369     {
370         return pme_lb->end;
371     }
372     else
373     {
374         return pme_lb->n;
375     }
376 }
377
378 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
379                                   gmx_int64_t step,
380                                   pme_load_balancing_t pme_lb)
381 {
382     char buf[STRLEN], sbuf[22];
383
384     sprintf(buf, "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
385             gmx_step_str(step, sbuf),
386             pmelblim_str[pme_lb->elimited],
387             pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
388     if (fp_err != NULL)
389     {
390         fprintf(fp_err, "\r%s\n", buf);
391     }
392     if (fp_log != NULL)
393     {
394         fprintf(fp_log, "%s\n", buf);
395     }
396 }
397
398 static void switch_to_stage1(pme_load_balancing_t pme_lb)
399 {
400     pme_lb->start = 0;
401     while (pme_lb->start+1 < pme_lb->n &&
402            (pme_lb->setup[pme_lb->start].count == 0 ||
403             pme_lb->setup[pme_lb->start].cycles >
404             pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC))
405     {
406         pme_lb->start++;
407     }
408     while (pme_lb->start > 0 && pme_lb->setup[pme_lb->start-1].cycles == 0)
409     {
410         pme_lb->start--;
411     }
412
413     pme_lb->end = pme_lb->n;
414     if (pme_lb->setup[pme_lb->end-1].count > 0 &&
415         pme_lb->setup[pme_lb->end-1].cycles >
416         pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
417     {
418         pme_lb->end--;
419     }
420
421     pme_lb->stage = 1;
422
423     /* Next we want to choose setup pme_lb->start, but as we will increase
424      * pme_ln->cur by one right after returning, we subtract 1 here.
425      */
426     pme_lb->cur = pme_lb->start - 1;
427 }
428
429 gmx_bool pme_load_balance(pme_load_balancing_t pme_lb,
430                           t_commrec           *cr,
431                           FILE                *fp_err,
432                           FILE                *fp_log,
433                           t_inputrec          *ir,
434                           t_state             *state,
435                           double               cycles,
436                           interaction_const_t *ic,
437                           nonbonded_verlet_t  *nbv,
438                           gmx_pme_t           *pmedata,
439                           gmx_int64_t          step)
440 {
441     gmx_bool     OK;
442     pme_setup_t *set;
443     double       cycles_fast;
444     char         buf[STRLEN], sbuf[22];
445     real         rtab;
446     gmx_bool     bUsesSimpleTables = TRUE;
447
448     if (pme_lb->stage == pme_lb->nstage)
449     {
450         return FALSE;
451     }
452
453     if (PAR(cr))
454     {
455         gmx_sumd(1, &cycles, cr);
456         cycles /= cr->nnodes;
457     }
458
459     set = &pme_lb->setup[pme_lb->cur];
460     set->count++;
461
462     rtab = ir->rlistlong + ir->tabext;
463
464     if (set->count % 2 == 1)
465     {
466         /* Skip the first cycle, because the first step after a switch
467          * is much slower due to allocation and/or caching effects.
468          */
469         return TRUE;
470     }
471
472     sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
473     print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
474
475     if (set->count <= 2)
476     {
477         set->cycles = cycles;
478     }
479     else
480     {
481         if (cycles*PME_LB_ACCEL_TOL < set->cycles &&
482             pme_lb->stage == pme_lb->nstage - 1)
483         {
484             /* The performance went up a lot (due to e.g. DD load balancing).
485              * Add a stage, keep the minima, but rescan all setups.
486              */
487             pme_lb->nstage++;
488
489             if (debug)
490             {
491                 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
492                         "Increased the number stages to %d"
493                         " and ignoring the previous performance\n",
494                         set->grid[XX], set->grid[YY], set->grid[ZZ],
495                         cycles*1e-6, set->cycles*1e-6, PME_LB_ACCEL_TOL,
496                         pme_lb->nstage);
497             }
498         }
499         set->cycles = min(set->cycles, cycles);
500     }
501
502     if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
503     {
504         pme_lb->fastest = pme_lb->cur;
505
506         if (DOMAINDECOMP(cr))
507         {
508             /* We found a new fastest setting, ensure that with subsequent
509              * shorter cut-off's the dynamic load balancing does not make
510              * the use of the current cut-off impossible. This solution is
511              * a trade-off, as the PME load balancing and DD domain size
512              * load balancing can interact in complex ways.
513              * With the Verlet kernels, DD load imbalance will usually be
514              * mainly due to bonded interaction imbalance, which will often
515              * quickly push the domain boundaries beyond the limit for the
516              * optimal, PME load balanced, cut-off. But it could be that
517              * better overal performance can be obtained with a slightly
518              * shorter cut-off and better DD load balancing.
519              */
520             change_dd_dlb_cutoff_limit(cr);
521         }
522     }
523     cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
524
525     /* Check in stage 0 if we should stop scanning grids.
526      * Stop when the time is more than SLOW_FAC longer than the fastest.
527      */
528     if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
529         cycles > pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
530     {
531         pme_lb->n = pme_lb->cur + 1;
532         /* Done with scanning, go to stage 1 */
533         switch_to_stage1(pme_lb);
534     }
535
536     if (pme_lb->stage == 0)
537     {
538         int gridsize_start;
539
540         gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
541
542         do
543         {
544             if (pme_lb->cur+1 < pme_lb->n)
545             {
546                 /* We had already generated the next setup */
547                 OK = TRUE;
548             }
549             else
550             {
551                 /* Find the next setup */
552                 OK = pme_loadbal_increase_cutoff(pme_lb, ir->pme_order, cr->dd);
553
554                 if (!OK)
555                 {
556                     pme_lb->elimited = epmelblimPMEGRID;
557                 }
558             }
559
560             if (OK && ir->ePBC != epbcNONE)
561             {
562                 OK = (sqr(pme_lb->setup[pme_lb->cur+1].rlistlong)
563                       <= max_cutoff2(ir->ePBC, state->box));
564                 if (!OK)
565                 {
566                     pme_lb->elimited = epmelblimBOX;
567                 }
568             }
569
570             if (OK)
571             {
572                 pme_lb->cur++;
573
574                 if (DOMAINDECOMP(cr))
575                 {
576                     OK = change_dd_cutoff(cr, state, ir,
577                                           pme_lb->setup[pme_lb->cur].rlistlong);
578                     if (!OK)
579                     {
580                         /* Failed: do not use this setup */
581                         pme_lb->cur--;
582                         pme_lb->elimited = epmelblimDD;
583                     }
584                 }
585             }
586             if (!OK)
587             {
588                 /* We hit the upper limit for the cut-off,
589                  * the setup should not go further than cur.
590                  */
591                 pme_lb->n = pme_lb->cur + 1;
592                 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
593                 /* Switch to the next stage */
594                 switch_to_stage1(pme_lb);
595             }
596         }
597         while (OK &&
598                !(pme_lb->setup[pme_lb->cur].grid[XX]*
599                  pme_lb->setup[pme_lb->cur].grid[YY]*
600                  pme_lb->setup[pme_lb->cur].grid[ZZ] <
601                  gridsize_start*PME_LB_GRID_SCALE_FAC
602                  &&
603                  pme_lb->setup[pme_lb->cur].grid_efficiency <
604                  pme_lb->setup[pme_lb->cur-1].grid_efficiency*PME_LB_GRID_EFFICIENCY_REL_FAC));
605     }
606
607     if (pme_lb->stage > 0 && pme_lb->end == 1)
608     {
609         pme_lb->cur   = 0;
610         pme_lb->stage = pme_lb->nstage;
611     }
612     else if (pme_lb->stage > 0 && pme_lb->end > 1)
613     {
614         /* If stage = nstage-1:
615          *   scan over all setups, rerunning only those setups
616          *   which are not much slower than the fastest
617          * else:
618          *   use the next setup
619          */
620         do
621         {
622             pme_lb->cur++;
623             if (pme_lb->cur == pme_lb->end)
624             {
625                 pme_lb->stage++;
626                 pme_lb->cur = pme_lb->start;
627             }
628         }
629         while (pme_lb->stage == pme_lb->nstage - 1 &&
630                pme_lb->setup[pme_lb->cur].count > 0 &&
631                pme_lb->setup[pme_lb->cur].cycles > cycles_fast*PME_LB_SLOW_FAC);
632
633         if (pme_lb->stage == pme_lb->nstage)
634         {
635             /* We are done optimizing, use the fastest setup we found */
636             pme_lb->cur = pme_lb->fastest;
637         }
638     }
639
640     if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
641     {
642         OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlistlong);
643         if (!OK)
644         {
645             /* Failsafe solution */
646             if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
647             {
648                 pme_lb->stage--;
649             }
650             pme_lb->fastest  = 0;
651             pme_lb->start    = 0;
652             pme_lb->end      = pme_lb->cur;
653             pme_lb->cur      = pme_lb->start;
654             pme_lb->elimited = epmelblimDD;
655             print_loadbal_limited(fp_err, fp_log, step, pme_lb);
656         }
657     }
658
659     /* Change the Coulomb cut-off and the PME grid */
660
661     set = &pme_lb->setup[pme_lb->cur];
662
663     ic->rcoulomb     = set->rcut_coulomb;
664     ic->rlist        = set->rlist;
665     ic->rlistlong    = set->rlistlong;
666     ir->nstcalclr    = set->nstcalclr;
667     ic->ewaldcoeff_q = set->ewaldcoeff_q;
668     /* TODO: centralize the code that sets the potentials shifts */
669     if (ic->coulomb_modifier == eintmodPOTSHIFT)
670     {
671         ic->sh_ewald = gmx_erfc(ic->ewaldcoeff_q*ic->rcoulomb);
672     }
673     if (EVDW_PME(ic->vdwtype))
674     {
675         /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
676         ic->rvdw            = set->rcut_coulomb;
677         ic->ewaldcoeff_lj   = set->ewaldcoeff_lj;
678         if (ic->vdw_modifier == eintmodPOTSHIFT)
679         {
680             real crc2;
681
682             ic->dispersion_shift.cpot = -pow(ic->rvdw, -6.0);
683             ic->repulsion_shift.cpot  = -pow(ic->rvdw, -12.0);
684             ic->sh_invrc6             = -ic->dispersion_shift.cpot;
685             crc2                      = sqr(ic->ewaldcoeff_lj*ic->rvdw);
686             ic->sh_lj_ewald           = (exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)*pow(ic->rvdw, -6.0);
687         }
688     }
689
690     bUsesSimpleTables = uses_simple_tables(ir->cutoff_scheme, nbv, 0);
691     if (pme_lb->cutoff_scheme == ecutsVERLET &&
692         nbv->grp[0].kernel_type == nbnxnk8x8x8_CUDA)
693     {
694         nbnxn_cuda_pme_loadbal_update_param(nbv->cu_nbv, ic);
695
696         /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
697          * also sharing texture references. To keep the code simple, we don't
698          * treat texture references as shared resources, but this means that
699          * the coulomb_tab texture ref will get updated by multiple threads.
700          * Hence, to ensure that the non-bonded kernels don't start before all
701          * texture binding operations are finished, we need to wait for all ranks
702          * to arrive here before continuing.
703          *
704          * Note that we could omit this barrier if GPUs are not shared (or
705          * texture objects are used), but as this is initialization code, there
706          * is not point in complicating things.
707          */
708 #ifdef GMX_THREAD_MPI
709         if (PAR(cr))
710         {
711             gmx_barrier(cr);
712         }
713 #endif  /* GMX_THREAD_MPI */
714     }
715
716     /* Usually we won't need the simple tables with GPUs.
717      * But we do with hybrid acceleration and with free energy.
718      * To avoid bugs, we always re-initialize the simple tables here.
719      */
720     init_interaction_const_tables(NULL, ic, bUsesSimpleTables, rtab);
721
722     if (cr->duty & DUTY_PME)
723     {
724         if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
725         {
726             /* Generate a new PME data structure,
727              * copying part of the old pointers.
728              */
729             gmx_pme_reinit(&set->pmedata,
730                            cr, pme_lb->setup[0].pmedata, ir,
731                            set->grid);
732         }
733         *pmedata = set->pmedata;
734     }
735     else
736     {
737         /* Tell our PME-only node to switch grid */
738         gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
739     }
740
741     if (debug)
742     {
743         print_grid(NULL, debug, "", "switched to", set, -1);
744     }
745
746     if (pme_lb->stage == pme_lb->nstage)
747     {
748         print_grid(fp_err, fp_log, "", "optimal", set, -1);
749     }
750
751     return TRUE;
752 }
753
754 void restart_pme_loadbal(pme_load_balancing_t pme_lb, int n)
755 {
756     pme_lb->nstage += n;
757 }
758
759 static int pme_grid_points(const pme_setup_t *setup)
760 {
761     return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
762 }
763
764 static real pme_loadbal_rlist(const pme_setup_t *setup)
765 {
766     /* With the group cut-off scheme we can have twin-range either
767      * for Coulomb or for VdW, so we need a check here.
768      * With the Verlet cut-off scheme rlist=rlistlong.
769      */
770     if (setup->rcut_coulomb > setup->rlist)
771     {
772         return setup->rlistlong;
773     }
774     else
775     {
776         return setup->rlist;
777     }
778 }
779
780 static void print_pme_loadbal_setting(FILE              *fplog,
781                                       char              *name,
782                                       const pme_setup_t *setup)
783 {
784     fprintf(fplog,
785             "   %-7s %6.3f nm %6.3f nm     %3d %3d %3d   %5.3f nm  %5.3f nm\n",
786             name,
787             setup->rcut_coulomb, pme_loadbal_rlist(setup),
788             setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
789             setup->spacing, 1/setup->ewaldcoeff_q);
790 }
791
792 static void print_pme_loadbal_settings(pme_load_balancing_t pme_lb,
793                                        t_commrec           *cr,
794                                        FILE                *fplog,
795                                        gmx_bool             bNonBondedOnGPU)
796 {
797     double pp_ratio, grid_ratio;
798
799     pp_ratio   = pow(pme_loadbal_rlist(&pme_lb->setup[pme_lb->cur])/pme_loadbal_rlist(&pme_lb->setup[0]), 3.0);
800     grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
801         (double)pme_grid_points(&pme_lb->setup[0]);
802
803     fprintf(fplog, "\n");
804     fprintf(fplog, "       P P   -   P M E   L O A D   B A L A N C I N G\n");
805     fprintf(fplog, "\n");
806     /* Here we only warn when the optimal setting is the last one */
807     if (pme_lb->elimited != epmelblimNO &&
808         pme_lb->cur == pme_loadbal_end(pme_lb)-1)
809     {
810         fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
811                 pmelblim_str[pme_lb->elimited]);
812         fprintf(fplog, "       you might not have reached a good load balance.\n");
813         if (pme_lb->elimited == epmelblimDD)
814         {
815             fprintf(fplog, "       Try different mdrun -dd settings or lower the -dds value.\n");
816         }
817         fprintf(fplog, "\n");
818     }
819     fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
820     fprintf(fplog, "           particle-particle                    PME\n");
821     fprintf(fplog, "            rcoulomb  rlist            grid      spacing   1/beta\n");
822     print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
823     print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
824     fprintf(fplog, " cost-ratio           %4.2f             %4.2f\n",
825             pp_ratio, grid_ratio);
826     fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
827
828     if (pp_ratio > 1.5 && !bNonBondedOnGPU)
829     {
830         md_print_warn(cr, fplog,
831                       "NOTE: PME load balancing increased the non-bonded workload by more than 50%%.\n"
832                       "      For better performance use (more) PME nodes (mdrun -npme),\n"
833                       "      or in case you are beyond the scaling limit, use less nodes in total.\n");
834     }
835     else
836     {
837         fprintf(fplog, "\n");
838     }
839 }
840
841 void pme_loadbal_done(pme_load_balancing_t pme_lb,
842                       t_commrec *cr, FILE *fplog,
843                       gmx_bool bNonBondedOnGPU)
844 {
845     if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
846     {
847         print_pme_loadbal_settings(pme_lb, cr, fplog, bNonBondedOnGPU);
848     }
849
850     /* TODO: Here we should free all pointers in pme_lb,
851      * but as it contains pme data structures,
852      * we need to first make pme.c free all data.
853      */
854 }