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