Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / mdlib / domdec_setup.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  * This file is part of Gromacs        Copyright (c) 1991-2008
5  * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
6  *
7  * This program is free software; you can redistribute it and/or
8  * modify it under the terms of the GNU General Public License
9  * as published by the Free Software Foundation; either version 2
10  * of the License, or (at your option) any later version.
11  *
12  * To help us fund GROMACS development, we humbly ask that you cite
13  * the research papers on the package. Check out http://www.gromacs.org
14  *
15  * And Hey:
16  * Gnomes, ROck Monsters And Chili Sauce
17  */
18
19 #ifdef HAVE_CONFIG_H
20 #include <config.h>
21 #endif
22
23 #include <stdio.h>
24 #include <assert.h>
25 #include "domdec.h"
26 #include "network.h"
27 #include "perf_est.h"
28 #include "physics.h"
29 #include "smalloc.h"
30 #include "typedefs.h"
31 #include "vec.h"
32 #include "names.h"
33
34 /* Margin for setting up the DD grid */
35 #define DD_GRID_MARGIN_PRES_SCALE 1.05
36
37 static int factorize(int n, int **fac, int **mfac)
38 {
39     int d, ndiv;
40
41     if (n <= 0)
42     {
43         gmx_fatal(FARGS, "Can only factorize positive integers.");
44     }
45
46     /* Decompose n in factors */
47     snew(*fac, n/2);
48     snew(*mfac, n/2);
49     d    = 2;
50     ndiv = 0;
51     while (n > 1)
52     {
53         while (n % d == 0)
54         {
55             if (ndiv == 0 || (*fac)[ndiv-1] != d)
56             {
57                 ndiv++;
58                 (*fac)[ndiv-1] = d;
59             }
60             (*mfac)[ndiv-1]++;
61             n /= d;
62         }
63         d++;
64     }
65
66     return ndiv;
67 }
68
69 static gmx_bool largest_divisor(int n)
70 {
71     int ndiv, *div, *mdiv, ldiv;
72
73     ndiv = factorize(n, &div, &mdiv);
74     ldiv = div[ndiv-1];
75     sfree(div);
76     sfree(mdiv);
77
78     return ldiv;
79 }
80
81 static int lcd(int n1, int n2)
82 {
83     int d, i;
84
85     d = 1;
86     for (i = 2; (i <= n1 && i <= n2); i++)
87     {
88         if (n1 % i == 0 && n2 % i == 0)
89         {
90             d = i;
91         }
92     }
93
94     return d;
95 }
96
97 static gmx_bool fits_pme_ratio(int nnodes, int npme, float ratio)
98 {
99     return ((double)npme/(double)nnodes > 0.95*ratio);
100 }
101
102 static gmx_bool fits_pp_pme_perf(int nnodes, int npme, float ratio)
103 {
104     int ndiv, *div, *mdiv, ldiv;
105     int npp_root3, npme_root2;
106
107     ndiv = factorize(nnodes-npme, &div, &mdiv);
108     ldiv = div[ndiv-1];
109     sfree(div);
110     sfree(mdiv);
111
112     npp_root3  = (int)(pow(nnodes-npme, 1.0/3.0) + 0.5);
113     npme_root2 = (int)(sqrt(npme) + 0.5);
114
115     /* The check below gives a reasonable division:
116      * factor 5 allowed at 5 or more PP nodes,
117      * factor 7 allowed at 49 or more PP nodes.
118      */
119     if (ldiv > 3 + npp_root3)
120     {
121         return FALSE;
122     }
123
124     /* Check if the number of PP and PME nodes have a reasonable sized
125      * denominator in common, such that we can use 2D PME decomposition
126      * when required (which requires nx_pp == nx_pme).
127      * The factor of 2 allows for a maximum ratio of 2^2=4
128      * between nx_pme and ny_pme.
129      */
130     if (lcd(nnodes-npme, npme)*2 < npme_root2)
131     {
132         return FALSE;
133     }
134
135     /* Does this division gives a reasonable PME load? */
136     return fits_pme_ratio(nnodes, npme, ratio);
137 }
138
139 static int guess_npme(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir, matrix box,
140                       int nnodes)
141 {
142     float      ratio;
143     int        npme, nkx, nky;
144     t_inputrec ir_try;
145
146     ratio = pme_load_estimate(mtop, ir, box);
147
148     if (fplog)
149     {
150         fprintf(fplog, "Guess for relative PME load: %.2f\n", ratio);
151     }
152
153     /* We assume the optimal node ratio is close to the load ratio.
154      * The communication load is neglected,
155      * but (hopefully) this will balance out between PP and PME.
156      */
157
158     if (!fits_pme_ratio(nnodes, nnodes/2, ratio))
159     {
160         /* We would need more than nnodes/2 PME only nodes,
161          * which is not possible. Since the PME load is very high,
162          * we will not loose much performance when all nodes do PME.
163          */
164
165         return 0;
166     }
167
168     /* First try to find npme as a factor of nnodes up to nnodes/3.
169      * We start with a minimum PME node fraction of 1/16
170      * and avoid ratios which lead to large prime factors in nnodes-npme.
171      */
172     npme = (nnodes + 15)/16;
173     while (npme <= nnodes/3)
174     {
175         if (nnodes % npme == 0)
176         {
177             /* Note that fits_perf might change the PME grid,
178              * in the current implementation it does not.
179              */
180             if (fits_pp_pme_perf(nnodes, npme, ratio))
181             {
182                 break;
183             }
184         }
185         npme++;
186     }
187     if (npme > nnodes/3)
188     {
189         /* Try any possible number for npme */
190         npme = 1;
191         while (npme <= nnodes/2)
192         {
193             /* Note that fits_perf may change the PME grid */
194             if (fits_pp_pme_perf(nnodes, npme, ratio))
195             {
196                 break;
197             }
198             npme++;
199         }
200     }
201     if (npme > nnodes/2)
202     {
203         gmx_fatal(FARGS, "Could not find an appropriate number of separate PME nodes. i.e. >= %5f*#nodes (%d) and <= #nodes/2 (%d) and reasonable performance wise (grid_x=%d, grid_y=%d).\n"
204                   "Use the -npme option of mdrun or change the number of processors or the PME grid dimensions, see the manual for details.",
205                   ratio, (int)(0.95*ratio*nnodes+0.5), nnodes/2, ir->nkx, ir->nky);
206         /* Keep the compiler happy */
207         npme = 0;
208     }
209     else
210     {
211         if (fplog)
212         {
213             fprintf(fplog,
214                     "Will use %d particle-particle and %d PME only nodes\n"
215                     "This is a guess, check the performance at the end of the log file\n",
216                     nnodes-npme, npme);
217         }
218         fprintf(stderr, "\n"
219                 "Will use %d particle-particle and %d PME only nodes\n"
220                 "This is a guess, check the performance at the end of the log file\n",
221                 nnodes-npme, npme);
222     }
223
224     return npme;
225 }
226
227 static int div_up(int n, int f)
228 {
229     return (n + f - 1)/f;
230 }
231
232 real comm_box_frac(ivec dd_nc, real cutoff, gmx_ddbox_t *ddbox)
233 {
234     int  i, j, k, npp;
235     rvec bt, nw;
236     real comm_vol;
237
238     for (i = 0; i < DIM; i++)
239     {
240         bt[i] = ddbox->box_size[i]*ddbox->skew_fac[i];
241         nw[i] = dd_nc[i]*cutoff/bt[i];
242     }
243
244     npp      = 1;
245     comm_vol = 0;
246     for (i = 0; i < DIM; i++)
247     {
248         if (dd_nc[i] > 1)
249         {
250             npp      *= dd_nc[i];
251             comm_vol += nw[i];
252             for (j = i+1; j < DIM; j++)
253             {
254                 if (dd_nc[j] > 1)
255                 {
256                     comm_vol += nw[i]*nw[j]*M_PI/4;
257                     for (k = j+1; k < DIM; k++)
258                     {
259                         if (dd_nc[k] > 1)
260                         {
261                             comm_vol += nw[i]*nw[j]*nw[k]*M_PI/6;
262                         }
263                     }
264                 }
265             }
266         }
267     }
268
269     return comm_vol;
270 }
271
272 static gmx_bool inhomogeneous_z(const t_inputrec *ir)
273 {
274     return ((EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD) &&
275             ir->ePBC == epbcXYZ && ir->ewald_geometry == eewg3DC);
276 }
277
278 /* Avoid integer overflows */
279 static float comm_pme_cost_vol(int npme, int a, int b, int c)
280 {
281     float comm_vol;
282
283     comm_vol  = npme - 1;
284     comm_vol *= npme;
285     comm_vol *= div_up(a, npme);
286     comm_vol *= div_up(b, npme);
287     comm_vol *= c;
288     return comm_vol;
289 }
290
291 static float comm_cost_est(real limit, real cutoff,
292                            matrix box, gmx_ddbox_t *ddbox,
293                            int natoms, t_inputrec *ir,
294                            float pbcdxr,
295                            int npme_tot, ivec nc)
296 {
297     ivec  npme = {1, 1, 1};
298     int   i, j, k, nk, overlap;
299     rvec  bt;
300     float comm_vol, comm_vol_xf, comm_pme, cost_pbcdx;
301     /* This is the cost of a pbc_dx call relative to the cost
302      * of communicating the coordinate and force of an atom.
303      * This will be machine dependent.
304      * These factors are for x86 with SMP or Infiniband.
305      */
306     float pbcdx_rect_fac = 0.1;
307     float pbcdx_tric_fac = 0.2;
308     float temp;
309
310     /* Check the DD algorithm restrictions */
311     if ((ir->ePBC == epbcXY && ir->nwall < 2 && nc[ZZ] > 1) ||
312         (ir->ePBC == epbcSCREW && (nc[XX] == 1 || nc[YY] > 1 || nc[ZZ] > 1)))
313     {
314         return -1;
315     }
316
317     if (inhomogeneous_z(ir) && nc[ZZ] > 1)
318     {
319         return -1;
320     }
321
322     assert(ddbox->npbcdim <= DIM);
323
324     /* Check if the triclinic requirements are met */
325     for (i = 0; i < DIM; i++)
326     {
327         for (j = i+1; j < ddbox->npbcdim; j++)
328         {
329             if (box[j][i] != 0 || ir->deform[j][i] != 0 ||
330                 (ir->epc != epcNO && ir->compress[j][i] != 0))
331             {
332                 if (nc[j] > 1 && nc[i] == 1)
333                 {
334                     return -1;
335                 }
336             }
337         }
338     }
339
340     for (i = 0; i < DIM; i++)
341     {
342         bt[i] = ddbox->box_size[i]*ddbox->skew_fac[i];
343
344         /* Without PBC there are no cell size limits with 2 cells */
345         if (!(i >= ddbox->npbcdim && nc[i] <= 2) && bt[i] < nc[i]*limit)
346         {
347             return -1;
348         }
349     }
350
351     if (npme_tot > 1)
352     {
353         /* The following choices should match those
354          * in init_domain_decomposition in domdec.c.
355          */
356         if (nc[XX] == 1 && nc[YY] > 1)
357         {
358             npme[XX] = 1;
359             npme[YY] = npme_tot;
360         }
361         else if (nc[YY] == 1)
362         {
363             npme[XX] = npme_tot;
364             npme[YY] = 1;
365         }
366         else
367         {
368             /* Will we use 1D or 2D PME decomposition? */
369             npme[XX] = (npme_tot % nc[XX] == 0) ? nc[XX] : npme_tot;
370             npme[YY] = npme_tot/npme[XX];
371         }
372     }
373
374     /* When two dimensions are (nearly) equal, use more cells
375      * for the smallest index, so the decomposition does not
376      * depend sensitively on the rounding of the box elements.
377      */
378     for (i = 0; i < DIM; i++)
379     {
380         for (j = i+1; j < DIM; j++)
381         {
382             /* Check if the box size is nearly identical,
383              * in that case we prefer nx > ny  and ny > nz.
384              */
385             if (fabs(bt[j] - bt[i]) < 0.01*bt[i] && nc[j] > nc[i])
386             {
387                 /* The XX/YY check is a bit compact. If nc[YY]==npme[YY]
388                  * this means the swapped nc has nc[XX]==npme[XX],
389                  * and we can also swap X and Y for PME.
390                  */
391                 /* Check if dimension i and j are equivalent for PME.
392                  * For x/y: if nc[YY]!=npme[YY], we can not swap x/y
393                  * For y/z: we can not have PME decomposition in z
394                  */
395                 if (npme_tot <= 1 ||
396                     !((i == XX && j == YY && nc[YY] != npme[YY]) ||
397                       (i == YY && j == ZZ && npme[YY] > 1)))
398                 {
399                     return -1;
400                 }
401             }
402         }
403     }
404
405     /* This function determines only half of the communication cost.
406      * All PP, PME and PP-PME communication is symmetric
407      * and the "back"-communication cost is identical to the forward cost.
408      */
409
410     comm_vol = comm_box_frac(nc, cutoff, ddbox);
411
412     comm_pme = 0;
413     for (i = 0; i < 2; i++)
414     {
415         /* Determine the largest volume for PME x/f redistribution */
416         if (nc[i] % npme[i] != 0)
417         {
418             if (nc[i] > npme[i])
419             {
420                 comm_vol_xf = (npme[i] == 2 ? 1.0/3.0 : 0.5);
421             }
422             else
423             {
424                 comm_vol_xf = 1.0 - lcd(nc[i], npme[i])/(double)npme[i];
425             }
426             comm_pme += 3*natoms*comm_vol_xf;
427         }
428
429         /* Grid overlap communication */
430         if (npme[i] > 1)
431         {
432             nk        = (i == 0 ? ir->nkx : ir->nky);
433             overlap   = (nk % npme[i] == 0 ? ir->pme_order-1 : ir->pme_order);
434             temp      = npme[i];
435             temp     *= overlap;
436             temp     *= ir->nkx;
437             temp     *= ir->nky;
438             temp     *= ir->nkz;
439             temp     /= nk;
440             comm_pme += temp;
441 /* Old line comm_pme += npme[i]*overlap*ir->nkx*ir->nky*ir->nkz/nk; */
442         }
443     }
444
445     /* PME FFT communication volume.
446      * This only takes the communication into account and not imbalance
447      * in the calculation. But the imbalance in communication and calculation
448      * are similar and therefore these formulas also prefer load balance
449      * in the FFT and pme_solve calculation.
450      */
451     comm_pme += comm_pme_cost_vol(npme[YY], ir->nky, ir->nkz, ir->nkx);
452     comm_pme += comm_pme_cost_vol(npme[XX], ir->nkx, ir->nky, ir->nkz);
453
454     /* Add cost of pbc_dx for bondeds */
455     cost_pbcdx = 0;
456     if ((nc[XX] == 1 || nc[YY] == 1) || (nc[ZZ] == 1 && ir->ePBC != epbcXY))
457     {
458         if ((ddbox->tric_dir[XX] && nc[XX] == 1) ||
459             (ddbox->tric_dir[YY] && nc[YY] == 1))
460         {
461             cost_pbcdx = pbcdxr*pbcdx_tric_fac;
462         }
463         else
464         {
465             cost_pbcdx = pbcdxr*pbcdx_rect_fac;
466         }
467     }
468
469     if (debug)
470     {
471         fprintf(debug,
472                 "nc %2d %2d %2d %2d %2d vol pp %6.4f pbcdx %6.4f pme %9.3e tot %9.3e\n",
473                 nc[XX], nc[YY], nc[ZZ], npme[XX], npme[YY],
474                 comm_vol, cost_pbcdx, comm_pme,
475                 3*natoms*(comm_vol + cost_pbcdx) + comm_pme);
476     }
477
478     return 3*natoms*(comm_vol + cost_pbcdx) + comm_pme;
479 }
480
481 static void assign_factors(gmx_domdec_t *dd,
482                            real limit, real cutoff,
483                            matrix box, gmx_ddbox_t *ddbox,
484                            int natoms, t_inputrec *ir,
485                            float pbcdxr, int npme,
486                            int ndiv, int *div, int *mdiv, ivec ir_try, ivec opt)
487 {
488     int   x, y, z, i;
489     float ce;
490
491     if (ndiv == 0)
492     {
493         ce = comm_cost_est(limit, cutoff, box, ddbox,
494                            natoms, ir, pbcdxr, npme, ir_try);
495         if (ce >= 0 && (opt[XX] == 0 ||
496                         ce < comm_cost_est(limit, cutoff, box, ddbox,
497                                            natoms, ir, pbcdxr,
498                                            npme, opt)))
499         {
500             copy_ivec(ir_try, opt);
501         }
502
503         return;
504     }
505
506     for (x = mdiv[0]; x >= 0; x--)
507     {
508         for (i = 0; i < x; i++)
509         {
510             ir_try[XX] *= div[0];
511         }
512         for (y = mdiv[0]-x; y >= 0; y--)
513         {
514             for (i = 0; i < y; i++)
515             {
516                 ir_try[YY] *= div[0];
517             }
518             for (i = 0; i < mdiv[0]-x-y; i++)
519             {
520                 ir_try[ZZ] *= div[0];
521             }
522
523             /* recurse */
524             assign_factors(dd, limit, cutoff, box, ddbox, natoms, ir, pbcdxr, npme,
525                            ndiv-1, div+1, mdiv+1, ir_try, opt);
526
527             for (i = 0; i < mdiv[0]-x-y; i++)
528             {
529                 ir_try[ZZ] /= div[0];
530             }
531             for (i = 0; i < y; i++)
532             {
533                 ir_try[YY] /= div[0];
534             }
535         }
536         for (i = 0; i < x; i++)
537         {
538             ir_try[XX] /= div[0];
539         }
540     }
541 }
542
543 static real optimize_ncells(FILE *fplog,
544                             int nnodes_tot, int npme_only,
545                             gmx_bool bDynLoadBal, real dlb_scale,
546                             gmx_mtop_t *mtop, matrix box, gmx_ddbox_t *ddbox,
547                             t_inputrec *ir,
548                             gmx_domdec_t *dd,
549                             real cellsize_limit, real cutoff,
550                             gmx_bool bInterCGBondeds,
551                             ivec nc)
552 {
553     int      npp, npme, ndiv, *div, *mdiv, d, nmax;
554     gmx_bool bExcl_pbcdx;
555     float    pbcdxr;
556     real     limit;
557     ivec     itry;
558
559     limit  = cellsize_limit;
560
561     dd->nc[XX] = 1;
562     dd->nc[YY] = 1;
563     dd->nc[ZZ] = 1;
564
565     npp = nnodes_tot - npme_only;
566     if (EEL_PME(ir->coulombtype))
567     {
568         npme = (npme_only > 0 ? npme_only : npp);
569     }
570     else
571     {
572         npme = 0;
573     }
574
575     if (bInterCGBondeds)
576     {
577         /* For Ewald exclusions pbc_dx is not called */
578         bExcl_pbcdx =
579             (IR_EXCL_FORCES(*ir) && !EEL_FULL(ir->coulombtype));
580         pbcdxr = (double)n_bonded_dx(mtop, bExcl_pbcdx)/(double)mtop->natoms;
581     }
582     else
583     {
584         /* Every molecule is a single charge group: no pbc required */
585         pbcdxr = 0;
586     }
587     /* Add a margin for DLB and/or pressure scaling */
588     if (bDynLoadBal)
589     {
590         if (dlb_scale >= 1.0)
591         {
592             gmx_fatal(FARGS, "The value for option -dds should be smaller than 1");
593         }
594         if (fplog)
595         {
596             fprintf(fplog, "Scaling the initial minimum size with 1/%g (option -dds) = %g\n", dlb_scale, 1/dlb_scale);
597         }
598         limit /= dlb_scale;
599     }
600     else if (ir->epc != epcNO)
601     {
602         if (fplog)
603         {
604             fprintf(fplog, "To account for pressure scaling, scaling the initial minimum size with %g\n", DD_GRID_MARGIN_PRES_SCALE);
605             limit *= DD_GRID_MARGIN_PRES_SCALE;
606         }
607     }
608
609     if (fplog)
610     {
611         fprintf(fplog, "Optimizing the DD grid for %d cells with a minimum initial size of %.3f nm\n", npp, limit);
612
613         if (inhomogeneous_z(ir))
614         {
615             fprintf(fplog, "Ewald_geometry=%s: assuming inhomogeneous particle distribution in z, will not decompose in z.\n", eewg_names[ir->ewald_geometry]);
616         }
617
618         if (limit > 0)
619         {
620             fprintf(fplog, "The maximum allowed number of cells is:");
621             for (d = 0; d < DIM; d++)
622             {
623                 nmax = (int)(ddbox->box_size[d]*ddbox->skew_fac[d]/limit);
624                 if (d >= ddbox->npbcdim && nmax < 2)
625                 {
626                     nmax = 2;
627                 }
628                 if (d == ZZ && inhomogeneous_z(ir))
629                 {
630                     nmax = 1;
631                 }
632                 fprintf(fplog, " %c %d", 'X' + d, nmax);
633             }
634             fprintf(fplog, "\n");
635         }
636     }
637
638     if (debug)
639     {
640         fprintf(debug, "Average nr of pbc_dx calls per atom %.2f\n", pbcdxr);
641     }
642
643     /* Decompose npp in factors */
644     ndiv = factorize(npp, &div, &mdiv);
645
646     itry[XX] = 1;
647     itry[YY] = 1;
648     itry[ZZ] = 1;
649     clear_ivec(nc);
650     assign_factors(dd, limit, cutoff, box, ddbox, mtop->natoms, ir, pbcdxr,
651                    npme, ndiv, div, mdiv, itry, nc);
652
653     sfree(div);
654     sfree(mdiv);
655
656     return limit;
657 }
658
659 real dd_choose_grid(FILE *fplog,
660                     t_commrec *cr, gmx_domdec_t *dd, t_inputrec *ir,
661                     gmx_mtop_t *mtop, matrix box, gmx_ddbox_t *ddbox,
662                     gmx_bool bDynLoadBal, real dlb_scale,
663                     real cellsize_limit, real cutoff_dd,
664                     gmx_bool bInterCGBondeds)
665 {
666     gmx_large_int_t nnodes_div, ldiv;
667     real            limit;
668
669     if (MASTER(cr))
670     {
671         nnodes_div = cr->nnodes;
672         if (EEL_PME(ir->coulombtype))
673         {
674             if (cr->npmenodes > 0)
675             {
676                 if (cr->nnodes <= 2)
677                 {
678                     gmx_fatal(FARGS,
679                               "Can not have separate PME nodes with 2 or less nodes");
680                 }
681                 if (cr->npmenodes >= cr->nnodes)
682                 {
683                     gmx_fatal(FARGS,
684                               "Can not have %d separate PME nodes with just %d total nodes",
685                               cr->npmenodes, cr->nnodes);
686                 }
687
688                 /* If the user purposely selected the number of PME nodes,
689                  * only check for large primes in the PP node count.
690                  */
691                 nnodes_div -= cr->npmenodes;
692             }
693         }
694         else
695         {
696             cr->npmenodes = 0;
697         }
698
699         if (cr->nnodes > 12)
700         {
701             ldiv = largest_divisor(nnodes_div);
702             /* Check if the largest divisor is more than nnodes^2/3 */
703             if (ldiv*ldiv*ldiv > nnodes_div*nnodes_div)
704             {
705                 gmx_fatal(FARGS, "The number of nodes you selected (%d) contains a large prime factor %d. In most cases this will lead to bad performance. Choose a number with smaller prime factors or set the decomposition (option -dd) manually.",
706                           nnodes_div, ldiv);
707             }
708         }
709
710         if (EEL_PME(ir->coulombtype))
711         {
712             if (cr->npmenodes < 0)
713             {
714                 /* Use PME nodes when the number of nodes is more than 16 */
715                 if (cr->nnodes <= 18)
716                 {
717                     cr->npmenodes = 0;
718                     if (fplog)
719                     {
720                         fprintf(fplog, "Using %d separate PME nodes, as there are too few total\n nodes for efficient splitting\n", cr->npmenodes);
721                     }
722                 }
723                 else
724                 {
725                     cr->npmenodes = guess_npme(fplog, mtop, ir, box, cr->nnodes);
726                     if (fplog)
727                     {
728                         fprintf(fplog, "Using %d separate PME nodes, as guessed by mdrun\n", cr->npmenodes);
729                     }
730                 }
731             }
732             else
733             {
734                 if (fplog)
735                 {
736                     fprintf(fplog, "Using %d separate PME nodes, per user request\n", cr->npmenodes);
737                 }
738             }
739         }
740
741         limit = optimize_ncells(fplog, cr->nnodes, cr->npmenodes,
742                                 bDynLoadBal, dlb_scale,
743                                 mtop, box, ddbox, ir, dd,
744                                 cellsize_limit, cutoff_dd,
745                                 bInterCGBondeds,
746                                 dd->nc);
747     }
748     else
749     {
750         limit = 0;
751     }
752     /* Communicate the information set by the master to all nodes */
753     gmx_bcast(sizeof(dd->nc), dd->nc, cr);
754     if (EEL_PME(ir->coulombtype))
755     {
756         gmx_bcast(sizeof(ir->nkx), &ir->nkx, cr);
757         gmx_bcast(sizeof(ir->nky), &ir->nky, cr);
758         gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);
759     }
760     else
761     {
762         cr->npmenodes = 0;
763     }
764
765     return limit;
766 }