Merge release-4-6 into master
[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, epmelblimNR
83 };
84
85 const char *pmelblim_str[epmelblimNR] =
86 { "no", "box size", "domain decompostion" };
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 {
209     pme_setup_t *set;
210     real         fac, sp;
211     real         tmpr_coulomb, tmpr_vdw;
212     int          d;
213
214     /* Try to add a new setup with next larger cut-off to the list */
215     pme_lb->n++;
216     srenew(pme_lb->setup, pme_lb->n);
217     set          = &pme_lb->setup[pme_lb->n-1];
218     set->pmedata = NULL;
219
220     fac = 1;
221     do
222     {
223         fac *= 1.01;
224         clear_ivec(set->grid);
225         sp = calc_grid(NULL, pme_lb->box_start,
226                        fac*pme_lb->setup[pme_lb->cur].spacing,
227                        &set->grid[XX],
228                        &set->grid[YY],
229                        &set->grid[ZZ]);
230
231         /* In parallel we can't have grids smaller than 2*pme_order,
232          * and we would anyhow not gain much speed at these grid sizes.
233          */
234         for (d = 0; d < DIM; d++)
235         {
236             if (set->grid[d] <= 2*pme_order)
237             {
238                 pme_lb->n--;
239
240                 return FALSE;
241             }
242         }
243     }
244     while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing);
245
246     set->rcut_coulomb = pme_lb->cut_spacing*sp;
247
248     if (pme_lb->cutoff_scheme == ecutsVERLET)
249     {
250         set->rlist        = set->rcut_coulomb + pme_lb->rbuf_coulomb;
251         /* We dont use LR lists with Verlet, but this avoids if-statements in further checks */
252         set->rlistlong    = set->rlist;
253     }
254     else
255     {
256         tmpr_coulomb          = set->rcut_coulomb + pme_lb->rbuf_coulomb;
257         tmpr_vdw              = pme_lb->rcut_vdw + pme_lb->rbuf_vdw;
258         set->rlist            = min(tmpr_coulomb, tmpr_vdw);
259         set->rlistlong        = max(tmpr_coulomb, tmpr_vdw);
260
261         /* Set the long-range update frequency */
262         if (set->rlist == set->rlistlong)
263         {
264             /* No long-range interactions if the short-/long-range cutoffs are identical */
265             set->nstcalclr = 0;
266         }
267         else if (pme_lb->nstcalclr_start == 0 || pme_lb->nstcalclr_start == 1)
268         {
269             /* We were not doing long-range before, but now we are since rlist!=rlistlong */
270             set->nstcalclr = 1;
271         }
272         else
273         {
274             /* We were already doing long-range interactions from the start */
275             if (pme_lb->rcut_vdw > pme_lb->rcut_coulomb_start)
276             {
277                 /* We were originally doing long-range VdW-only interactions.
278                  * If rvdw is still longer than rcoulomb we keep the original nstcalclr,
279                  * but if the coulomb cutoff has become longer we should update the long-range
280                  * part every step.
281                  */
282                 set->nstcalclr = (tmpr_vdw > tmpr_coulomb) ? pme_lb->nstcalclr_start : 1;
283             }
284             else
285             {
286                 /* We were not doing any long-range interaction from the start,
287                  * since it is not possible to do twin-range coulomb for the PME interaction.
288                  */
289                 set->nstcalclr = 1;
290             }
291         }
292     }
293
294     set->spacing      = sp;
295     /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
296     set->grid_efficiency = 1;
297     for (d = 0; d < DIM; d++)
298     {
299         set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
300     }
301     /* The Ewald coefficient is inversly proportional to the cut-off */
302     set->ewaldcoeff =
303         pme_lb->setup[0].ewaldcoeff*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
304
305     set->count   = 0;
306     set->cycles  = 0;
307
308     if (debug)
309     {
310         fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
311                 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
312     }
313     return TRUE;
314 }
315
316 static void print_grid(FILE *fp_err, FILE *fp_log,
317                        const char *pre,
318                        const char *desc,
319                        const pme_setup_t *set,
320                        double cycles)
321 {
322     char buf[STRLEN], buft[STRLEN];
323
324     if (cycles >= 0)
325     {
326         sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
327     }
328     else
329     {
330         buft[0] = '\0';
331     }
332     sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
333             pre,
334             desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
335             buft);
336     if (fp_err != NULL)
337     {
338         fprintf(fp_err, "\r%s\n", buf);
339     }
340     if (fp_log != NULL)
341     {
342         fprintf(fp_log, "%s\n", buf);
343     }
344 }
345
346 static int pme_loadbal_end(pme_load_balancing_t pme_lb)
347 {
348     /* In the initial stage only n is set; end is not set yet */
349     if (pme_lb->end > 0)
350     {
351         return pme_lb->end;
352     }
353     else
354     {
355         return pme_lb->n;
356     }
357 }
358
359 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
360                                   gmx_large_int_t step,
361                                   pme_load_balancing_t pme_lb)
362 {
363     char buf[STRLEN], sbuf[22];
364
365     sprintf(buf, "step %4s: the %s limited the PME load balancing to a coulomb cut-off of %.3f",
366             gmx_step_str(step, sbuf),
367             pmelblim_str[pme_lb->elimited],
368             pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
369     if (fp_err != NULL)
370     {
371         fprintf(fp_err, "\r%s\n", buf);
372     }
373     if (fp_log != NULL)
374     {
375         fprintf(fp_log, "%s\n", buf);
376     }
377 }
378
379 static void switch_to_stage1(pme_load_balancing_t pme_lb)
380 {
381     pme_lb->start = 0;
382     while (pme_lb->start+1 < pme_lb->n &&
383            (pme_lb->setup[pme_lb->start].count == 0 ||
384             pme_lb->setup[pme_lb->start].cycles >
385             pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC))
386     {
387         pme_lb->start++;
388     }
389     while (pme_lb->start > 0 && pme_lb->setup[pme_lb->start-1].cycles == 0)
390     {
391         pme_lb->start--;
392     }
393
394     pme_lb->end = pme_lb->n;
395     if (pme_lb->setup[pme_lb->end-1].count > 0 &&
396         pme_lb->setup[pme_lb->end-1].cycles >
397         pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
398     {
399         pme_lb->end--;
400     }
401
402     pme_lb->stage = 1;
403
404     /* Next we want to choose setup pme_lb->start, but as we will increase
405      * pme_ln->cur by one right after returning, we subtract 1 here.
406      */
407     pme_lb->cur = pme_lb->start - 1;
408 }
409
410 gmx_bool pme_load_balance(pme_load_balancing_t pme_lb,
411                           t_commrec           *cr,
412                           FILE                *fp_err,
413                           FILE                *fp_log,
414                           t_inputrec          *ir,
415                           t_state             *state,
416                           double               cycles,
417                           interaction_const_t *ic,
418                           nonbonded_verlet_t  *nbv,
419                           gmx_pme_t           *pmedata,
420                           gmx_large_int_t      step)
421 {
422     gmx_bool     OK;
423     pme_setup_t *set;
424     double       cycles_fast;
425     char         buf[STRLEN], sbuf[22];
426     real         rtab;
427     gmx_bool     bUsesSimpleTables = TRUE;
428
429     if (pme_lb->stage == pme_lb->nstage)
430     {
431         return FALSE;
432     }
433
434     if (PAR(cr))
435     {
436         gmx_sumd(1, &cycles, cr);
437         cycles /= cr->nnodes;
438     }
439
440     set = &pme_lb->setup[pme_lb->cur];
441     set->count++;
442
443     rtab = ir->rlistlong + ir->tabext;
444
445     if (set->count % 2 == 1)
446     {
447         /* Skip the first cycle, because the first step after a switch
448          * is much slower due to allocation and/or caching effects.
449          */
450         return TRUE;
451     }
452
453     sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
454     print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
455
456     if (set->count <= 2)
457     {
458         set->cycles = cycles;
459     }
460     else
461     {
462         if (cycles*PME_LB_ACCEL_TOL < set->cycles &&
463             pme_lb->stage == pme_lb->nstage - 1)
464         {
465             /* The performance went up a lot (due to e.g. DD load balancing).
466              * Add a stage, keep the minima, but rescan all setups.
467              */
468             pme_lb->nstage++;
469
470             if (debug)
471             {
472                 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
473                         "Increased the number stages to %d"
474                         " and ignoring the previous performance\n",
475                         set->grid[XX], set->grid[YY], set->grid[ZZ],
476                         cycles*1e-6, set->cycles*1e-6, PME_LB_ACCEL_TOL,
477                         pme_lb->nstage);
478             }
479         }
480         set->cycles = min(set->cycles, cycles);
481     }
482
483     if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
484     {
485         pme_lb->fastest = pme_lb->cur;
486
487         if (DOMAINDECOMP(cr))
488         {
489             /* We found a new fastest setting, ensure that with subsequent
490              * shorter cut-off's the dynamic load balancing does not make
491              * the use of the current cut-off impossible. This solution is
492              * a trade-off, as the PME load balancing and DD domain size
493              * load balancing can interact in complex ways.
494              * With the Verlet kernels, DD load imbalance will usually be
495              * mainly due to bonded interaction imbalance, which will often
496              * quickly push the domain boundaries beyond the limit for the
497              * optimal, PME load balanced, cut-off. But it could be that
498              * better overal performance can be obtained with a slightly
499              * shorter cut-off and better DD load balancing.
500              */
501             change_dd_dlb_cutoff_limit(cr);
502         }
503     }
504     cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
505
506     /* Check in stage 0 if we should stop scanning grids.
507      * Stop when the time is more than SLOW_FAC longer than the fastest.
508      */
509     if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
510         cycles > pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
511     {
512         pme_lb->n = pme_lb->cur + 1;
513         /* Done with scanning, go to stage 1 */
514         switch_to_stage1(pme_lb);
515     }
516
517     if (pme_lb->stage == 0)
518     {
519         int gridsize_start;
520
521         gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
522
523         do
524         {
525             if (pme_lb->cur+1 < pme_lb->n)
526             {
527                 /* We had already generated the next setup */
528                 OK = TRUE;
529             }
530             else
531             {
532                 /* Find the next setup */
533                 OK = pme_loadbal_increase_cutoff(pme_lb, ir->pme_order);
534             }
535
536             if (OK && ir->ePBC != epbcNONE)
537             {
538                 OK = (sqr(pme_lb->setup[pme_lb->cur+1].rlistlong)
539                       <= max_cutoff2(ir->ePBC, state->box));
540                 if (!OK)
541                 {
542                     pme_lb->elimited = epmelblimBOX;
543                 }
544             }
545
546             if (OK)
547             {
548                 pme_lb->cur++;
549
550                 if (DOMAINDECOMP(cr))
551                 {
552                     OK = change_dd_cutoff(cr, state, ir,
553                                           pme_lb->setup[pme_lb->cur].rlistlong);
554                     if (!OK)
555                     {
556                         /* Failed: do not use this setup */
557                         pme_lb->cur--;
558                         pme_lb->elimited = epmelblimDD;
559                     }
560                 }
561             }
562             if (!OK)
563             {
564                 /* We hit the upper limit for the cut-off,
565                  * the setup should not go further than cur.
566                  */
567                 pme_lb->n = pme_lb->cur + 1;
568                 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
569                 /* Switch to the next stage */
570                 switch_to_stage1(pme_lb);
571             }
572         }
573         while (OK &&
574                !(pme_lb->setup[pme_lb->cur].grid[XX]*
575                  pme_lb->setup[pme_lb->cur].grid[YY]*
576                  pme_lb->setup[pme_lb->cur].grid[ZZ] <
577                  gridsize_start*PME_LB_GRID_SCALE_FAC
578                  &&
579                  pme_lb->setup[pme_lb->cur].grid_efficiency <
580                  pme_lb->setup[pme_lb->cur-1].grid_efficiency*PME_LB_GRID_EFFICIENCY_REL_FAC));
581     }
582
583     if (pme_lb->stage > 0 && pme_lb->end == 1)
584     {
585         pme_lb->cur   = 0;
586         pme_lb->stage = pme_lb->nstage;
587     }
588     else if (pme_lb->stage > 0 && pme_lb->end > 1)
589     {
590         /* If stage = nstage-1:
591          *   scan over all setups, rerunning only those setups
592          *   which are not much slower than the fastest
593          * else:
594          *   use the next setup
595          */
596         do
597         {
598             pme_lb->cur++;
599             if (pme_lb->cur == pme_lb->end)
600             {
601                 pme_lb->stage++;
602                 pme_lb->cur = pme_lb->start;
603             }
604         }
605         while (pme_lb->stage == pme_lb->nstage - 1 &&
606                pme_lb->setup[pme_lb->cur].count > 0 &&
607                pme_lb->setup[pme_lb->cur].cycles > cycles_fast*PME_LB_SLOW_FAC);
608
609         if (pme_lb->stage == pme_lb->nstage)
610         {
611             /* We are done optimizing, use the fastest setup we found */
612             pme_lb->cur = pme_lb->fastest;
613         }
614     }
615
616     if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
617     {
618         OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlistlong);
619         if (!OK)
620         {
621             /* Failsafe solution */
622             if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
623             {
624                 pme_lb->stage--;
625             }
626             pme_lb->fastest  = 0;
627             pme_lb->start    = 0;
628             pme_lb->end      = pme_lb->cur;
629             pme_lb->cur      = pme_lb->start;
630             pme_lb->elimited = epmelblimDD;
631             print_loadbal_limited(fp_err, fp_log, step, pme_lb);
632         }
633     }
634
635     /* Change the Coulomb cut-off and the PME grid */
636
637     set = &pme_lb->setup[pme_lb->cur];
638
639     ic->rcoulomb   = set->rcut_coulomb;
640     ic->rlist      = set->rlist;
641     ic->rlistlong  = set->rlistlong;
642     ir->nstcalclr  = set->nstcalclr;
643     ic->ewaldcoeff = set->ewaldcoeff;
644
645     bUsesSimpleTables = uses_simple_tables(ir->cutoff_scheme, nbv, 0);
646     if (pme_lb->cutoff_scheme == ecutsVERLET &&
647         nbv->grp[0].kernel_type == nbnxnk8x8x8_CUDA)
648     {
649         nbnxn_cuda_pme_loadbal_update_param(nbv->cu_nbv, ic);
650     }
651     else
652     {
653         init_interaction_const_tables(NULL, ic, bUsesSimpleTables,
654                                       rtab);
655     }
656
657     if (pme_lb->cutoff_scheme == ecutsVERLET && nbv->ngrp > 1)
658     {
659         init_interaction_const_tables(NULL, ic, bUsesSimpleTables,
660                                       rtab);
661     }
662
663     if (cr->duty & DUTY_PME)
664     {
665         if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
666         {
667             /* Generate a new PME data structure,
668              * copying part of the old pointers.
669              */
670             gmx_pme_reinit(&set->pmedata,
671                            cr, pme_lb->setup[0].pmedata, ir,
672                            set->grid);
673         }
674         *pmedata = set->pmedata;
675     }
676     else
677     {
678         /* Tell our PME-only node to switch grid */
679         gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff);
680     }
681
682     if (debug)
683     {
684         print_grid(NULL, debug, "", "switched to", set, -1);
685     }
686
687     if (pme_lb->stage == pme_lb->nstage)
688     {
689         print_grid(fp_err, fp_log, "", "optimal", set, -1);
690     }
691
692     return TRUE;
693 }
694
695 void restart_pme_loadbal(pme_load_balancing_t pme_lb, int n)
696 {
697     pme_lb->nstage += n;
698 }
699
700 static int pme_grid_points(const pme_setup_t *setup)
701 {
702     return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
703 }
704
705 static real pme_loadbal_rlist(const pme_setup_t *setup)
706 {
707     /* With the group cut-off scheme we can have twin-range either
708      * for Coulomb or for VdW, so we need a check here.
709      * With the Verlet cut-off scheme rlist=rlistlong.
710      */
711     if (setup->rcut_coulomb > setup->rlist)
712     {
713         return setup->rlistlong;
714     }
715     else
716     {
717         return setup->rlist;
718     }
719 }
720
721 static void print_pme_loadbal_setting(FILE              *fplog,
722                                       char              *name,
723                                       const pme_setup_t *setup)
724 {
725     fprintf(fplog,
726             "   %-7s %6.3f nm %6.3f nm     %3d %3d %3d   %5.3f nm  %5.3f nm\n",
727             name,
728             setup->rcut_coulomb, pme_loadbal_rlist(setup),
729             setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
730             setup->spacing, 1/setup->ewaldcoeff);
731 }
732
733 static void print_pme_loadbal_settings(pme_load_balancing_t pme_lb,
734                                        t_commrec           *cr,
735                                        FILE                *fplog,
736                                        gmx_bool             bNonBondedOnGPU)
737 {
738     double pp_ratio, grid_ratio;
739
740     pp_ratio   = pow(pme_loadbal_rlist(&pme_lb->setup[pme_lb->cur])/pme_loadbal_rlist(&pme_lb->setup[0]), 3.0);
741     grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
742         (double)pme_grid_points(&pme_lb->setup[0]);
743
744     fprintf(fplog, "\n");
745     fprintf(fplog, "       P P   -   P M E   L O A D   B A L A N C I N G\n");
746     fprintf(fplog, "\n");
747     /* Here we only warn when the optimal setting is the last one */
748     if (pme_lb->elimited != epmelblimNO &&
749         pme_lb->cur == pme_loadbal_end(pme_lb)-1)
750     {
751         fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
752                 pmelblim_str[pme_lb->elimited]);
753         fprintf(fplog, "       you might not have reached a good load balance.\n");
754         if (pme_lb->elimited == epmelblimDD)
755         {
756             fprintf(fplog, "       Try different mdrun -dd settings or lower the -dds value.\n");
757         }
758         fprintf(fplog, "\n");
759     }
760     fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
761     fprintf(fplog, "           particle-particle                    PME\n");
762     fprintf(fplog, "            rcoulomb  rlist            grid      spacing   1/beta\n");
763     print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
764     print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
765     fprintf(fplog, " cost-ratio           %4.2f             %4.2f\n",
766             pp_ratio, grid_ratio);
767     fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
768
769     if (pp_ratio > 1.5 && !bNonBondedOnGPU)
770     {
771         md_print_warn(cr, fplog,
772                       "NOTE: PME load balancing increased the non-bonded workload by more than 50%%.\n"
773                       "      For better performance use (more) PME nodes (mdrun -npme),\n"
774                       "      or in case you are beyond the scaling limit, use less nodes in total.\n");
775     }
776     else
777     {
778         fprintf(fplog, "\n");
779     }
780 }
781
782 void pme_loadbal_done(pme_load_balancing_t pme_lb,
783                       t_commrec *cr, FILE *fplog,
784                       gmx_bool bNonBondedOnGPU)
785 {
786     if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
787     {
788         print_pme_loadbal_settings(pme_lb, cr, fplog, bNonBondedOnGPU);
789     }
790
791     /* TODO: Here we should free all pointers in pme_lb,
792      * but as it contains pme data structures,
793      * we need to first make pme.c free all data.
794      */
795 }