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