Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / programs / mdrun / pme_loadbal.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 4.6.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2011, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include "smalloc.h"
41 #include "network.h"
42 #include "calcgrid.h"
43 #include "pme.h"
44 #include "vec.h"
45 #include "domdec.h"
46 #include "nbnxn_cuda_data_mgmt.h"
47 #include "force.h"
48 #include "macros.h"
49 #include "md_logging.h"
50 #include "pme_loadbal.h"
51
52 /* Parameters and setting for one PP-PME setup */
53 typedef struct {
54     real      rcut_coulomb;    /* Coulomb cut-off                              */
55     real      rlist;           /* pair-list cut-off                            */
56     real      rlistlong;       /* LR pair-list cut-off                         */
57     int       nstcalclr;       /* frequency of evaluating long-range forces for group scheme */
58     real      spacing;         /* (largest) PME grid spacing                   */
59     ivec      grid;            /* the PME grid dimensions                      */
60     real      grid_efficiency; /* ineffiency factor for non-uniform grids <= 1 */
61     real      ewaldcoeff;      /* the Ewald coefficient                        */
62     gmx_pme_t pmedata;         /* the data structure used in the PME code      */
63
64     int       count;           /* number of times this setup has been timed    */
65     double    cycles;          /* the fastest time for this setup in cycles    */
66 } pme_setup_t;
67
68 /* In the initial scan, step by grids that are at least a factor 0.8 coarser */
69 #define PME_LB_GRID_SCALE_FAC  0.8
70 /* In the initial scan, try to skip grids with uneven x/y/z spacing,
71  * checking if the "efficiency" is more than 5% worse than the previous grid.
72  */
73 #define PME_LB_GRID_EFFICIENCY_REL_FAC  1.05
74 /* Rerun up till 12% slower setups than the fastest up till now */
75 #define PME_LB_SLOW_FAC  1.12
76 /* If setups get more than 2% faster, do another round to avoid
77  * choosing a slower setup due to acceleration or fluctuations.
78  */
79 #define PME_LB_ACCEL_TOL 1.02
80
81 enum {
82     epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimPMEGRID, epmelblimNR
83 };
84
85 const char *pmelblim_str[epmelblimNR] =
86 { "no", "box size", "domain decompostion", "PME grid restriction" };
87
88 struct pme_load_balancing {
89     int          nstage;             /* the current maximum number of stages */
90
91     real         cut_spacing;        /* the minimum cutoff / PME grid spacing ratio */
92     real         rcut_vdw;           /* Vdw cutoff (does not change) */
93     real         rcut_coulomb_start; /* Initial electrostatics cutoff */
94     int          nstcalclr_start;    /* Initial electrostatics cutoff */
95     real         rbuf_coulomb;       /* the pairlist buffer size */
96     real         rbuf_vdw;           /* the pairlist buffer size */
97     matrix       box_start;          /* the initial simulation box */
98     int          n;                  /* the count of setup as well as the allocation size */
99     pme_setup_t *setup;              /* the PME+cutoff setups */
100     int          cur;                /* the current setup */
101     int          fastest;            /* fastest setup up till now */
102     int          start;              /* start of setup range to consider in stage>0 */
103     int          end;                /* end   of setup range to consider in stage>0 */
104     int          elimited;           /* was the balancing limited, uses enum above */
105     int          cutoff_scheme;      /* Verlet or group cut-offs */
106
107     int          stage;              /* the current stage */
108 };
109
110 void pme_loadbal_init(pme_load_balancing_t *pme_lb_p,
111                       const t_inputrec *ir, matrix box,
112                       const interaction_const_t *ic,
113                       gmx_pme_t pmedata)
114 {
115     pme_load_balancing_t pme_lb;
116     real                 spm, sp;
117     int                  d;
118
119     snew(pme_lb, 1);
120
121     /* Any number of stages >= 2 is supported */
122     pme_lb->nstage   = 2;
123
124     pme_lb->cutoff_scheme = ir->cutoff_scheme;
125
126     if (pme_lb->cutoff_scheme == ecutsVERLET)
127     {
128         pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
129         pme_lb->rbuf_vdw     = pme_lb->rbuf_coulomb;
130     }
131     else
132     {
133         if (ic->rcoulomb > ic->rlist)
134         {
135             pme_lb->rbuf_coulomb = ic->rlistlong - ic->rcoulomb;
136         }
137         else
138         {
139             pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
140         }
141         if (ic->rvdw > ic->rlist)
142         {
143             pme_lb->rbuf_vdw = ic->rlistlong - ic->rvdw;
144         }
145         else
146         {
147             pme_lb->rbuf_vdw = ic->rlist - ic->rvdw;
148         }
149     }
150
151     copy_mat(box, pme_lb->box_start);
152     if (ir->ePBC == epbcXY && ir->nwall == 2)
153     {
154         svmul(ir->wall_ewald_zfac, pme_lb->box_start[ZZ], pme_lb->box_start[ZZ]);
155     }
156
157     pme_lb->n = 1;
158     snew(pme_lb->setup, pme_lb->n);
159
160     pme_lb->rcut_vdw              = ic->rvdw;
161     pme_lb->rcut_coulomb_start    = ir->rcoulomb;
162     pme_lb->nstcalclr_start       = ir->nstcalclr;
163
164     pme_lb->cur                   = 0;
165     pme_lb->setup[0].rcut_coulomb = ic->rcoulomb;
166     pme_lb->setup[0].rlist        = ic->rlist;
167     pme_lb->setup[0].rlistlong    = ic->rlistlong;
168     pme_lb->setup[0].nstcalclr    = ir->nstcalclr;
169     pme_lb->setup[0].grid[XX]     = ir->nkx;
170     pme_lb->setup[0].grid[YY]     = ir->nky;
171     pme_lb->setup[0].grid[ZZ]     = ir->nkz;
172     pme_lb->setup[0].ewaldcoeff   = ic->ewaldcoeff;
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 =
319         pme_lb->setup[0].ewaldcoeff*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
320
321     set->count   = 0;
322     set->cycles  = 0;
323
324     if (debug)
325     {
326         fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
327                 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
328     }
329     return TRUE;
330 }
331
332 static void print_grid(FILE *fp_err, FILE *fp_log,
333                        const char *pre,
334                        const char *desc,
335                        const pme_setup_t *set,
336                        double cycles)
337 {
338     char buf[STRLEN], buft[STRLEN];
339
340     if (cycles >= 0)
341     {
342         sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
343     }
344     else
345     {
346         buft[0] = '\0';
347     }
348     sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
349             pre,
350             desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
351             buft);
352     if (fp_err != NULL)
353     {
354         fprintf(fp_err, "\r%s\n", buf);
355     }
356     if (fp_log != NULL)
357     {
358         fprintf(fp_log, "%s\n", buf);
359     }
360 }
361
362 static int pme_loadbal_end(pme_load_balancing_t pme_lb)
363 {
364     /* In the initial stage only n is set; end is not set yet */
365     if (pme_lb->end > 0)
366     {
367         return pme_lb->end;
368     }
369     else
370     {
371         return pme_lb->n;
372     }
373 }
374
375 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
376                                   gmx_large_int_t step,
377                                   pme_load_balancing_t pme_lb)
378 {
379     char buf[STRLEN], sbuf[22];
380
381     sprintf(buf, "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
382             gmx_step_str(step, sbuf),
383             pmelblim_str[pme_lb->elimited],
384             pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
385     if (fp_err != NULL)
386     {
387         fprintf(fp_err, "\r%s\n", buf);
388     }
389     if (fp_log != NULL)
390     {
391         fprintf(fp_log, "%s\n", buf);
392     }
393 }
394
395 static void switch_to_stage1(pme_load_balancing_t pme_lb)
396 {
397     pme_lb->start = 0;
398     while (pme_lb->start+1 < pme_lb->n &&
399            (pme_lb->setup[pme_lb->start].count == 0 ||
400             pme_lb->setup[pme_lb->start].cycles >
401             pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC))
402     {
403         pme_lb->start++;
404     }
405     while (pme_lb->start > 0 && pme_lb->setup[pme_lb->start-1].cycles == 0)
406     {
407         pme_lb->start--;
408     }
409
410     pme_lb->end = pme_lb->n;
411     if (pme_lb->setup[pme_lb->end-1].count > 0 &&
412         pme_lb->setup[pme_lb->end-1].cycles >
413         pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
414     {
415         pme_lb->end--;
416     }
417
418     pme_lb->stage = 1;
419
420     /* Next we want to choose setup pme_lb->start, but as we will increase
421      * pme_ln->cur by one right after returning, we subtract 1 here.
422      */
423     pme_lb->cur = pme_lb->start - 1;
424 }
425
426 gmx_bool pme_load_balance(pme_load_balancing_t pme_lb,
427                           t_commrec           *cr,
428                           FILE                *fp_err,
429                           FILE                *fp_log,
430                           t_inputrec          *ir,
431                           t_state             *state,
432                           double               cycles,
433                           interaction_const_t *ic,
434                           nonbonded_verlet_t  *nbv,
435                           gmx_pme_t           *pmedata,
436                           gmx_large_int_t      step)
437 {
438     gmx_bool     OK;
439     pme_setup_t *set;
440     double       cycles_fast;
441     char         buf[STRLEN], sbuf[22];
442     real         rtab;
443     gmx_bool     bUsesSimpleTables = TRUE;
444
445     if (pme_lb->stage == pme_lb->nstage)
446     {
447         return FALSE;
448     }
449
450     if (PAR(cr))
451     {
452         gmx_sumd(1, &cycles, cr);
453         cycles /= cr->nnodes;
454     }
455
456     set = &pme_lb->setup[pme_lb->cur];
457     set->count++;
458
459     rtab = ir->rlistlong + ir->tabext;
460
461     if (set->count % 2 == 1)
462     {
463         /* Skip the first cycle, because the first step after a switch
464          * is much slower due to allocation and/or caching effects.
465          */
466         return TRUE;
467     }
468
469     sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
470     print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
471
472     if (set->count <= 2)
473     {
474         set->cycles = cycles;
475     }
476     else
477     {
478         if (cycles*PME_LB_ACCEL_TOL < set->cycles &&
479             pme_lb->stage == pme_lb->nstage - 1)
480         {
481             /* The performance went up a lot (due to e.g. DD load balancing).
482              * Add a stage, keep the minima, but rescan all setups.
483              */
484             pme_lb->nstage++;
485
486             if (debug)
487             {
488                 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
489                         "Increased the number stages to %d"
490                         " and ignoring the previous performance\n",
491                         set->grid[XX], set->grid[YY], set->grid[ZZ],
492                         cycles*1e-6, set->cycles*1e-6, PME_LB_ACCEL_TOL,
493                         pme_lb->nstage);
494             }
495         }
496         set->cycles = min(set->cycles, cycles);
497     }
498
499     if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
500     {
501         pme_lb->fastest = pme_lb->cur;
502
503         if (DOMAINDECOMP(cr))
504         {
505             /* We found a new fastest setting, ensure that with subsequent
506              * shorter cut-off's the dynamic load balancing does not make
507              * the use of the current cut-off impossible. This solution is
508              * a trade-off, as the PME load balancing and DD domain size
509              * load balancing can interact in complex ways.
510              * With the Verlet kernels, DD load imbalance will usually be
511              * mainly due to bonded interaction imbalance, which will often
512              * quickly push the domain boundaries beyond the limit for the
513              * optimal, PME load balanced, cut-off. But it could be that
514              * better overal performance can be obtained with a slightly
515              * shorter cut-off and better DD load balancing.
516              */
517             change_dd_dlb_cutoff_limit(cr);
518         }
519     }
520     cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
521
522     /* Check in stage 0 if we should stop scanning grids.
523      * Stop when the time is more than SLOW_FAC longer than the fastest.
524      */
525     if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
526         cycles > pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
527     {
528         pme_lb->n = pme_lb->cur + 1;
529         /* Done with scanning, go to stage 1 */
530         switch_to_stage1(pme_lb);
531     }
532
533     if (pme_lb->stage == 0)
534     {
535         int gridsize_start;
536
537         gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
538
539         do
540         {
541             if (pme_lb->cur+1 < pme_lb->n)
542             {
543                 /* We had already generated the next setup */
544                 OK = TRUE;
545             }
546             else
547             {
548                 /* Find the next setup */
549                 OK = pme_loadbal_increase_cutoff(pme_lb, ir->pme_order, cr->dd);
550
551                 if (!OK)
552                 {
553                     pme_lb->elimited = epmelblimPMEGRID;
554                 }
555             }
556
557             if (OK && ir->ePBC != epbcNONE)
558             {
559                 OK = (sqr(pme_lb->setup[pme_lb->cur+1].rlistlong)
560                       <= max_cutoff2(ir->ePBC, state->box));
561                 if (!OK)
562                 {
563                     pme_lb->elimited = epmelblimBOX;
564                 }
565             }
566
567             if (OK)
568             {
569                 pme_lb->cur++;
570
571                 if (DOMAINDECOMP(cr))
572                 {
573                     OK = change_dd_cutoff(cr, state, ir,
574                                           pme_lb->setup[pme_lb->cur].rlistlong);
575                     if (!OK)
576                     {
577                         /* Failed: do not use this setup */
578                         pme_lb->cur--;
579                         pme_lb->elimited = epmelblimDD;
580                     }
581                 }
582             }
583             if (!OK)
584             {
585                 /* We hit the upper limit for the cut-off,
586                  * the setup should not go further than cur.
587                  */
588                 pme_lb->n = pme_lb->cur + 1;
589                 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
590                 /* Switch to the next stage */
591                 switch_to_stage1(pme_lb);
592             }
593         }
594         while (OK &&
595                !(pme_lb->setup[pme_lb->cur].grid[XX]*
596                  pme_lb->setup[pme_lb->cur].grid[YY]*
597                  pme_lb->setup[pme_lb->cur].grid[ZZ] <
598                  gridsize_start*PME_LB_GRID_SCALE_FAC
599                  &&
600                  pme_lb->setup[pme_lb->cur].grid_efficiency <
601                  pme_lb->setup[pme_lb->cur-1].grid_efficiency*PME_LB_GRID_EFFICIENCY_REL_FAC));
602     }
603
604     if (pme_lb->stage > 0 && pme_lb->end == 1)
605     {
606         pme_lb->cur   = 0;
607         pme_lb->stage = pme_lb->nstage;
608     }
609     else if (pme_lb->stage > 0 && pme_lb->end > 1)
610     {
611         /* If stage = nstage-1:
612          *   scan over all setups, rerunning only those setups
613          *   which are not much slower than the fastest
614          * else:
615          *   use the next setup
616          */
617         do
618         {
619             pme_lb->cur++;
620             if (pme_lb->cur == pme_lb->end)
621             {
622                 pme_lb->stage++;
623                 pme_lb->cur = pme_lb->start;
624             }
625         }
626         while (pme_lb->stage == pme_lb->nstage - 1 &&
627                pme_lb->setup[pme_lb->cur].count > 0 &&
628                pme_lb->setup[pme_lb->cur].cycles > cycles_fast*PME_LB_SLOW_FAC);
629
630         if (pme_lb->stage == pme_lb->nstage)
631         {
632             /* We are done optimizing, use the fastest setup we found */
633             pme_lb->cur = pme_lb->fastest;
634         }
635     }
636
637     if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
638     {
639         OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlistlong);
640         if (!OK)
641         {
642             /* Failsafe solution */
643             if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
644             {
645                 pme_lb->stage--;
646             }
647             pme_lb->fastest  = 0;
648             pme_lb->start    = 0;
649             pme_lb->end      = pme_lb->cur;
650             pme_lb->cur      = pme_lb->start;
651             pme_lb->elimited = epmelblimDD;
652             print_loadbal_limited(fp_err, fp_log, step, pme_lb);
653         }
654     }
655
656     /* Change the Coulomb cut-off and the PME grid */
657
658     set = &pme_lb->setup[pme_lb->cur];
659
660     ic->rcoulomb   = set->rcut_coulomb;
661     ic->rlist      = set->rlist;
662     ic->rlistlong  = set->rlistlong;
663     ir->nstcalclr  = set->nstcalclr;
664     ic->ewaldcoeff = set->ewaldcoeff;
665
666     bUsesSimpleTables = uses_simple_tables(ir->cutoff_scheme, nbv, 0);
667     if (pme_lb->cutoff_scheme == ecutsVERLET &&
668         nbv->grp[0].kernel_type == nbnxnk8x8x8_CUDA)
669     {
670         nbnxn_cuda_pme_loadbal_update_param(nbv->cu_nbv, ic);
671     }
672     else
673     {
674         init_interaction_const_tables(NULL, ic, bUsesSimpleTables,
675                                       rtab);
676     }
677
678     if (pme_lb->cutoff_scheme == ecutsVERLET && nbv->ngrp > 1)
679     {
680         init_interaction_const_tables(NULL, ic, bUsesSimpleTables,
681                                       rtab);
682     }
683
684     if (cr->duty & DUTY_PME)
685     {
686         if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
687         {
688             /* Generate a new PME data structure,
689              * copying part of the old pointers.
690              */
691             gmx_pme_reinit(&set->pmedata,
692                            cr, pme_lb->setup[0].pmedata, ir,
693                            set->grid);
694         }
695         *pmedata = set->pmedata;
696     }
697     else
698     {
699         /* Tell our PME-only node to switch grid */
700         gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff);
701     }
702
703     if (debug)
704     {
705         print_grid(NULL, debug, "", "switched to", set, -1);
706     }
707
708     if (pme_lb->stage == pme_lb->nstage)
709     {
710         print_grid(fp_err, fp_log, "", "optimal", set, -1);
711     }
712
713     return TRUE;
714 }
715
716 void restart_pme_loadbal(pme_load_balancing_t pme_lb, int n)
717 {
718     pme_lb->nstage += n;
719 }
720
721 static int pme_grid_points(const pme_setup_t *setup)
722 {
723     return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
724 }
725
726 static real pme_loadbal_rlist(const pme_setup_t *setup)
727 {
728     /* With the group cut-off scheme we can have twin-range either
729      * for Coulomb or for VdW, so we need a check here.
730      * With the Verlet cut-off scheme rlist=rlistlong.
731      */
732     if (setup->rcut_coulomb > setup->rlist)
733     {
734         return setup->rlistlong;
735     }
736     else
737     {
738         return setup->rlist;
739     }
740 }
741
742 static void print_pme_loadbal_setting(FILE              *fplog,
743                                       char              *name,
744                                       const pme_setup_t *setup)
745 {
746     fprintf(fplog,
747             "   %-7s %6.3f nm %6.3f nm     %3d %3d %3d   %5.3f nm  %5.3f nm\n",
748             name,
749             setup->rcut_coulomb, pme_loadbal_rlist(setup),
750             setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
751             setup->spacing, 1/setup->ewaldcoeff);
752 }
753
754 static void print_pme_loadbal_settings(pme_load_balancing_t pme_lb,
755                                        t_commrec           *cr,
756                                        FILE                *fplog,
757                                        gmx_bool             bNonBondedOnGPU)
758 {
759     double pp_ratio, grid_ratio;
760
761     pp_ratio   = pow(pme_loadbal_rlist(&pme_lb->setup[pme_lb->cur])/pme_loadbal_rlist(&pme_lb->setup[0]), 3.0);
762     grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
763         (double)pme_grid_points(&pme_lb->setup[0]);
764
765     fprintf(fplog, "\n");
766     fprintf(fplog, "       P P   -   P M E   L O A D   B A L A N C I N G\n");
767     fprintf(fplog, "\n");
768     /* Here we only warn when the optimal setting is the last one */
769     if (pme_lb->elimited != epmelblimNO &&
770         pme_lb->cur == pme_loadbal_end(pme_lb)-1)
771     {
772         fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
773                 pmelblim_str[pme_lb->elimited]);
774         fprintf(fplog, "       you might not have reached a good load balance.\n");
775         if (pme_lb->elimited == epmelblimDD)
776         {
777             fprintf(fplog, "       Try different mdrun -dd settings or lower the -dds value.\n");
778         }
779         fprintf(fplog, "\n");
780     }
781     fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
782     fprintf(fplog, "           particle-particle                    PME\n");
783     fprintf(fplog, "            rcoulomb  rlist            grid      spacing   1/beta\n");
784     print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
785     print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
786     fprintf(fplog, " cost-ratio           %4.2f             %4.2f\n",
787             pp_ratio, grid_ratio);
788     fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
789
790     if (pp_ratio > 1.5 && !bNonBondedOnGPU)
791     {
792         md_print_warn(cr, fplog,
793                       "NOTE: PME load balancing increased the non-bonded workload by more than 50%%.\n"
794                       "      For better performance use (more) PME nodes (mdrun -npme),\n"
795                       "      or in case you are beyond the scaling limit, use less nodes in total.\n");
796     }
797     else
798     {
799         fprintf(fplog, "\n");
800     }
801 }
802
803 void pme_loadbal_done(pme_load_balancing_t pme_lb,
804                       t_commrec *cr, FILE *fplog,
805                       gmx_bool bNonBondedOnGPU)
806 {
807     if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
808     {
809         print_pme_loadbal_settings(pme_lb, cr, fplog, bNonBondedOnGPU);
810     }
811
812     /* TODO: Here we should free all pointers in pme_lb,
813      * but as it contains pme data structures,
814      * we need to first make pme.c free all data.
815      */
816 }