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