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