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