2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 * \brief This file defines functions used by the domdec module
39 * in its initial setup phase.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
47 #include "domdec_setup.h"
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_struct.h"
56 #include "gromacs/domdec/options.h"
57 #include "gromacs/ewald/pme.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/math/utilities.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/mdlib/perf_est.h"
62 #include "gromacs/mdtypes/commrec.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/topology/topology.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/logger.h"
69 #include "gromacs/utility/stringutil.h"
72 #include "domdec_internal.h"
75 // TODO remove this when moving domdec into gmx namespace
76 using gmx::DomdecOptions;
78 /*! \brief Margin for setting up the DD grid */
79 #define DD_GRID_MARGIN_PRES_SCALE 1.05
81 /*! \brief Factorize \p n.
83 * \param[in] n Value to factorize
84 * \param[out] fac Vector of factors (to be allocated in this function)
85 * \param[out] mfac Vector with the number of times each factor repeats in the factorization (to be allocated in this function)
87 static void factorize(int n, std::vector<int>* fac, std::vector<int>* mfac)
91 gmx_fatal(FARGS, "Can only factorize positive integers.");
94 /* Decompose n in factors */
102 if (fac->empty() || fac->back() != d)
117 /*! \brief Find largest divisor of \p n smaller than \p n*/
118 static int largest_divisor(int n)
120 std::vector<int> div;
121 std::vector<int> mdiv;
122 factorize(n, &div, &mdiv);
127 /*! \brief Compute largest common divisor of \p n1 and \b n2 */
128 static int lcd(int n1, int n2)
133 for (i = 2; (i <= n1 && i <= n2); i++)
135 if (n1 % i == 0 && n2 % i == 0)
144 /*! \brief Returns TRUE when there are enough PME ranks for the ratio */
145 static gmx_bool fits_pme_ratio(int nrank_tot, int nrank_pme, float ratio)
147 return (static_cast<double>(nrank_pme) / static_cast<double>(nrank_tot) > 0.95 * ratio);
150 /*! \brief Returns TRUE when npme out of ntot ranks doing PME is expected to give reasonable performance */
151 static gmx_bool fits_pp_pme_perf(int ntot, int npme, float ratio)
153 std::vector<int> div;
154 std::vector<int> mdiv;
155 factorize(ntot - npme, &div, &mdiv);
157 int npp_root3 = gmx::roundToInt(std::cbrt(ntot - npme));
158 int npme_root2 = gmx::roundToInt(std::sqrt(static_cast<double>(npme)));
160 /* The check below gives a reasonable division:
161 * factor 5 allowed at 5 or more PP ranks,
162 * factor 7 allowed at 49 or more PP ranks.
164 if (div.back() > 3 + npp_root3)
169 /* Check if the number of PP and PME ranks have a reasonable sized
170 * denominator in common, such that we can use 2D PME decomposition
171 * when required (which requires nx_pp == nx_pme).
172 * The factor of 2 allows for a maximum ratio of 2^2=4
173 * between nx_pme and ny_pme.
175 if (lcd(ntot - npme, npme) * 2 < npme_root2)
180 /* Does this division gives a reasonable PME load? */
181 return fits_pme_ratio(ntot, npme, ratio);
184 /*! \brief Make a guess for the number of PME ranks to use. */
185 static int guess_npme(const gmx::MDLogger& mdlog,
186 const gmx_mtop_t& mtop,
187 const t_inputrec& ir,
194 ratio = pme_load_estimate(mtop, ir, box);
196 GMX_LOG(mdlog.info).appendTextFormatted("Guess for relative PME load: %.2f", ratio);
198 /* We assume the optimal rank ratio is close to the load ratio.
199 * The communication load is neglected,
200 * but (hopefully) this will balance out between PP and PME.
203 if (!fits_pme_ratio(nrank_tot, nrank_tot / 2, ratio))
205 /* We would need more than nrank_tot/2 PME only nodes,
206 * which is not possible. Since the PME load is very high,
207 * we will not loose much performance when all ranks do PME.
213 /* First try to find npme as a factor of nrank_tot up to nrank_tot/3.
214 * We start with a minimum PME node fraction of 1/16
215 * and avoid ratios which lead to large prime factors in nnodes-npme.
217 npme = (nrank_tot + 15) / 16;
218 while (npme <= nrank_tot / 3)
220 if (nrank_tot % npme == 0)
222 /* Note that fits_perf might change the PME grid,
223 * in the current implementation it does not.
225 if (fits_pp_pme_perf(nrank_tot, npme, ratio))
232 if (npme > nrank_tot / 3)
234 /* Try any possible number for npme */
236 while (npme <= nrank_tot / 2)
238 /* Note that fits_perf may change the PME grid */
239 if (fits_pp_pme_perf(nrank_tot, npme, ratio))
246 if (npme > nrank_tot / 2)
249 "Could not find an appropriate number of separate PME ranks. i.e. >= %5f*#ranks "
250 "(%d) and <= #ranks/2 (%d) and reasonable performance wise (grid_x=%d, "
252 "Use the -npme option of mdrun or change the number of ranks or the PME grid "
253 "dimensions, see the manual for details.",
254 ratio, gmx::roundToInt(0.95 * ratio * nrank_tot), nrank_tot / 2, ir.nkx, ir.nky);
259 .appendTextFormatted(
260 "Will use %d particle-particle and %d PME only ranks\n"
261 "This is a guess, check the performance at the end of the log file",
262 nrank_tot - npme, npme);
268 /*! \brief Return \p n divided by \p f rounded up to the next integer. */
269 static int div_up(int n, int f)
271 return (n + f - 1) / f;
274 real comm_box_frac(const gmx::IVec& dd_nc, real cutoff, const gmx_ddbox_t& ddbox)
280 for (i = 0; i < DIM; i++)
282 real bt = ddbox.box_size[i] * ddbox.skew_fac[i];
283 nw[i] = dd_nc[i] * cutoff / bt;
287 for (i = 0; i < DIM; i++)
292 for (j = i + 1; j < DIM; j++)
296 comm_vol += nw[i] * nw[j] * M_PI / 4;
297 for (k = j + 1; k < DIM; k++)
301 comm_vol += nw[i] * nw[j] * nw[k] * M_PI / 6;
312 /*! \brief Return whether the DD inhomogeneous in the z direction */
313 static gmx_bool inhomogeneous_z(const t_inputrec& ir)
315 return ((EEL_PME(ir.coulombtype) || ir.coulombtype == eelEWALD) && ir.ePBC == epbcXYZ
316 && ir.ewald_geometry == eewg3DC);
319 /*! \brief Estimate cost of PME FFT communication
321 * This only takes the communication into account and not imbalance
322 * in the calculation. But the imbalance in communication and calculation
323 * are similar and therefore these formulas also prefer load balance
324 * in the FFT and pme_solve calculation.
326 static float comm_pme_cost_vol(int npme, int a, int b, int c)
328 /* We use a float here, since an integer might overflow */
333 comm_vol *= div_up(a, npme);
334 comm_vol *= div_up(b, npme);
340 /*! \brief Estimate cost of communication for a possible domain decomposition. */
341 static float comm_cost_est(real limit,
344 const gmx_ddbox_t& ddbox,
346 const t_inputrec& ir,
351 gmx::IVec npme = { 1, 1, 1 };
352 int i, j, nk, overlap;
354 float comm_vol, comm_vol_xf, comm_pme, cost_pbcdx;
355 /* This is the cost of a pbc_dx call relative to the cost
356 * of communicating the coordinate and force of an atom.
357 * This will be machine dependent.
358 * These factors are for x86 with SMP or Infiniband.
360 float pbcdx_rect_fac = 0.1;
361 float pbcdx_tric_fac = 0.2;
364 /* Check the DD algorithm restrictions */
365 if ((ir.ePBC == epbcXY && ir.nwall < 2 && nc[ZZ] > 1)
366 || (ir.ePBC == epbcSCREW && (nc[XX] == 1 || nc[YY] > 1 || nc[ZZ] > 1)))
371 if (inhomogeneous_z(ir) && nc[ZZ] > 1)
376 assert(ddbox.npbcdim <= DIM);
378 /* Check if the triclinic requirements are met */
379 for (i = 0; i < DIM; i++)
381 for (j = i + 1; j < ddbox.npbcdim; j++)
383 if (box[j][i] != 0 || ir.deform[j][i] != 0 || (ir.epc != epcNO && ir.compress[j][i] != 0))
385 if (nc[j] > 1 && nc[i] == 1)
393 for (i = 0; i < DIM; i++)
395 bt[i] = ddbox.box_size[i] * ddbox.skew_fac[i];
397 /* Without PBC and with 2 cells, there are no lower limits on the cell size */
398 if (!(i >= ddbox.npbcdim && nc[i] <= 2) && bt[i] < nc[i] * limit)
402 /* With PBC, check if the cut-off fits in nc[i]-1 cells */
403 if (i < ddbox.npbcdim && nc[i] > 1 && (nc[i] - 1) * bt[i] < nc[i] * cutoff)
411 /* The following choices should match those
412 * in init_domain_decomposition in domdec.c.
414 if (nc[XX] == 1 && nc[YY] > 1)
419 else if (nc[YY] == 1)
426 /* Will we use 1D or 2D PME decomposition? */
427 npme[XX] = (npme_tot % nc[XX] == 0) ? nc[XX] : npme_tot;
428 npme[YY] = npme_tot / npme[XX];
432 if (EEL_PME(ir.coulombtype) || EVDW_PME(ir.vdwtype))
434 /* Check the PME grid restrictions.
435 * Currently these can only be invalid here with too few grid lines
436 * along the x dimension per rank doing PME.
438 int npme_x = (npme_tot > 1 ? npme[XX] : nc[XX]);
440 /* Currently we don't have the OpenMP thread count available here.
441 * But with threads we have only tighter restrictions and it's
442 * probably better anyhow to avoid settings where we need to reduce
443 * grid lines over multiple ranks, as the thread check will do.
445 bool useThreads = true;
446 bool errorsAreFatal = false;
447 if (!gmx_pme_check_restrictions(ir.pme_order, ir.nkx, ir.nky, ir.nkz, npme_x, useThreads,
454 /* When two dimensions are (nearly) equal, use more cells
455 * for the smallest index, so the decomposition does not
456 * depend sensitively on the rounding of the box elements.
458 for (i = 0; i < DIM; i++)
460 for (j = i + 1; j < DIM; j++)
462 /* Check if the box size is nearly identical,
463 * in that case we prefer nx > ny and ny > nz.
465 if (std::fabs(bt[j] - bt[i]) < 0.01 * bt[i] && nc[j] > nc[i])
467 /* The XX/YY check is a bit compact. If nc[YY]==npme[YY]
468 * this means the swapped nc has nc[XX]==npme[XX],
469 * and we can also swap X and Y for PME.
471 /* Check if dimension i and j are equivalent for PME.
472 * For x/y: if nc[YY]!=npme[YY], we can not swap x/y
473 * For y/z: we can not have PME decomposition in z
476 || !((i == XX && j == YY && nc[YY] != npme[YY]) || (i == YY && j == ZZ && npme[YY] > 1)))
484 /* This function determines only half of the communication cost.
485 * All PP, PME and PP-PME communication is symmetric
486 * and the "back"-communication cost is identical to the forward cost.
489 comm_vol = comm_box_frac(nc, cutoff, ddbox);
492 for (i = 0; i < 2; i++)
494 /* Determine the largest volume for PME x/f redistribution */
495 if (nc[i] % npme[i] != 0)
499 comm_vol_xf = (npme[i] == 2 ? 1.0 / 3.0 : 0.5);
503 comm_vol_xf = 1.0 - lcd(nc[i], npme[i]) / static_cast<double>(npme[i]);
505 comm_pme += 3 * natoms * comm_vol_xf;
508 /* Grid overlap communication */
511 nk = (i == 0 ? ir.nkx : ir.nky);
512 overlap = (nk % npme[i] == 0 ? ir.pme_order - 1 : ir.pme_order);
520 /* Old line comm_pme += npme[i]*overlap*ir.nkx*ir.nky*ir.nkz/nk; */
524 comm_pme += comm_pme_cost_vol(npme[YY], ir.nky, ir.nkz, ir.nkx);
525 comm_pme += comm_pme_cost_vol(npme[XX], ir.nkx, ir.nky, ir.nkz);
527 /* Add cost of pbc_dx for bondeds */
529 if ((nc[XX] == 1 || nc[YY] == 1) || (nc[ZZ] == 1 && ir.ePBC != epbcXY))
531 if ((ddbox.tric_dir[XX] && nc[XX] == 1) || (ddbox.tric_dir[YY] && nc[YY] == 1))
533 cost_pbcdx = pbcdxr * pbcdx_tric_fac;
537 cost_pbcdx = pbcdxr * pbcdx_rect_fac;
543 fprintf(debug, "nc %2d %2d %2d %2d %2d vol pp %6.4f pbcdx %6.4f pme %9.3e tot %9.3e\n",
544 nc[XX], nc[YY], nc[ZZ], npme[XX], npme[YY], comm_vol, cost_pbcdx,
545 comm_pme / (3 * natoms), comm_vol + cost_pbcdx + comm_pme / (3 * natoms));
548 return 3 * natoms * (comm_vol + cost_pbcdx) + comm_pme;
551 /*! \brief Assign penalty factors to possible domain decompositions,
552 * based on the estimated communication costs. */
553 static void assign_factors(const real limit,
554 const bool request1D,
557 const gmx_ddbox_t& ddbox,
559 const t_inputrec& ir,
570 gmx::IVec& ir_try = *irTryPtr;
574 const int maxDimensionSize = std::max(ir_try[XX], std::max(ir_try[YY], ir_try[ZZ]));
575 const int productOfDimensionSizes = ir_try[XX] * ir_try[YY] * ir_try[ZZ];
576 const bool decompositionHasOneDimension = (maxDimensionSize == productOfDimensionSizes);
577 if (request1D && !decompositionHasOneDimension)
579 /* We requested 1D DD, but got multiple dimensions */
583 ce = comm_cost_est(limit, cutoff, box, ddbox, natoms, ir, pbcdxr, npme, ir_try);
586 || ce < comm_cost_est(limit, cutoff, box, ddbox, natoms, ir, pbcdxr, npme, *opt)))
594 for (x = mdiv[0]; x >= 0; x--)
596 for (i = 0; i < x; i++)
598 ir_try[XX] *= div[0];
600 for (y = mdiv[0] - x; y >= 0; y--)
602 for (i = 0; i < y; i++)
604 ir_try[YY] *= div[0];
606 for (i = 0; i < mdiv[0] - x - y; i++)
608 ir_try[ZZ] *= div[0];
612 assign_factors(limit, request1D, cutoff, box, ddbox, natoms, ir, pbcdxr, npme, ndiv - 1,
613 div + 1, mdiv + 1, irTryPtr, opt);
615 for (i = 0; i < mdiv[0] - x - y; i++)
617 ir_try[ZZ] /= div[0];
619 for (i = 0; i < y; i++)
621 ir_try[YY] /= div[0];
624 for (i = 0; i < x; i++)
626 ir_try[XX] /= div[0];
631 /*! \brief Determine the optimal distribution of DD cells for the
632 * simulation system and number of MPI ranks
634 * \returns The optimal grid cell choice. The latter will contain all
635 * zeros if no valid cell choice exists. */
636 static gmx::IVec optimizeDDCells(const gmx::MDLogger& mdlog,
637 const int numRanksRequested,
638 const int numPmeOnlyRanks,
639 const real cellSizeLimit,
640 const bool request1DAnd1Pulse,
641 const gmx_mtop_t& mtop,
643 const gmx_ddbox_t& ddbox,
644 const t_inputrec& ir,
645 const DDSystemInfo& systemInfo)
649 const int numPPRanks = numRanksRequested - numPmeOnlyRanks;
652 .appendTextFormatted(
653 "Optimizing the DD grid for %d cells with a minimum initial size of %.3f nm",
654 numPPRanks, cellSizeLimit);
655 if (inhomogeneous_z(ir))
658 .appendTextFormatted(
659 "Ewald_geometry=%s: assuming inhomogeneous particle distribution in z, "
660 "will not decompose in z.",
661 eewg_names[ir.ewald_geometry]);
665 // For cost estimates, we need the number of ranks doing PME work,
666 // which is the number of PP ranks when not using separate
668 const int numRanksDoingPmeWork =
669 (EEL_PME(ir.coulombtype) ? ((numPmeOnlyRanks > 0) ? numPmeOnlyRanks : numPPRanks) : 0);
671 if (systemInfo.haveInterDomainBondeds)
673 /* If we can skip PBC for distance calculations in plain-C bondeds,
674 * we can save some time (e.g. 3D DD with pbc=xyz).
675 * Here we ignore SIMD bondeds as they always do (fast) PBC.
677 count_bonded_distances(mtop, ir, &pbcdxr, nullptr);
678 pbcdxr /= static_cast<double>(mtop.natoms);
682 /* Every molecule is a single charge group: no pbc required */
686 if (cellSizeLimit > 0)
688 std::string maximumCells = "The maximum allowed number of cells is:";
689 for (int d = 0; d < DIM; d++)
691 int nmax = static_cast<int>(ddbox.box_size[d] * ddbox.skew_fac[d] / cellSizeLimit);
692 if (d >= ddbox.npbcdim && nmax < 2)
696 if (d == ZZ && inhomogeneous_z(ir))
700 maximumCells += gmx::formatString(" %c %d", 'X' + d, nmax);
702 GMX_LOG(mdlog.info).appendText(maximumCells);
707 fprintf(debug, "Average nr of pbc_dx calls per atom %.2f\n", pbcdxr);
710 /* Decompose numPPRanks in factors */
711 std::vector<int> div;
712 std::vector<int> mdiv;
713 factorize(numPPRanks, &div, &mdiv);
715 gmx::IVec itry = { 1, 1, 1 };
716 gmx::IVec numDomains = { 0, 0, 0 };
717 assign_factors(cellSizeLimit, request1DAnd1Pulse, systemInfo.cutoff, box, ddbox, mtop.natoms, ir,
718 pbcdxr, numRanksDoingPmeWork, div.size(), div.data(), mdiv.data(), &itry, &numDomains);
723 real getDDGridSetupCellSizeLimit(const gmx::MDLogger& mdlog,
724 const bool request1DAnd1Pulse,
725 const bool bDynLoadBal,
726 const real dlb_scale,
727 const t_inputrec& ir,
728 const real systemInfoCellSizeLimit)
730 real cellSizeLimit = systemInfoCellSizeLimit;
731 if (request1DAnd1Pulse)
733 cellSizeLimit = std::max(cellSizeLimit, ir.rlist);
736 /* Add a margin for DLB and/or pressure scaling */
739 if (dlb_scale >= 1.0)
741 gmx_fatal(FARGS, "The value for option -dds should be smaller than 1");
744 .appendTextFormatted(
745 "Scaling the initial minimum size with 1/%g (option -dds) = %g", dlb_scale,
747 cellSizeLimit /= dlb_scale;
749 else if (ir.epc != epcNO)
752 .appendTextFormatted(
753 "To account for pressure scaling, scaling the initial minimum size with %g",
754 DD_GRID_MARGIN_PRES_SCALE);
755 cellSizeLimit *= DD_GRID_MARGIN_PRES_SCALE;
758 return cellSizeLimit;
760 void checkForValidRankCountRequests(const int numRanksRequested, const bool usingPme, const int numPmeRanksRequested)
762 int numPPRanksRequested = numRanksRequested;
763 if (usingPme && numPmeRanksRequested > 0)
765 numPPRanksRequested -= numPmeRanksRequested;
766 if (numPmeRanksRequested > numPPRanksRequested)
769 "Cannot have %d separate PME ranks with only %d PP ranks, choose fewer or no "
770 "separate PME ranks",
771 numPmeRanksRequested, numPPRanksRequested);
775 // Once the rank count is large enough, it becomes worth
776 // suggesting improvements to the user.
777 const int minPPRankCountToCheckForLargePrimeFactors = 13;
778 if (numPPRanksRequested >= minPPRankCountToCheckForLargePrimeFactors)
780 const int largestDivisor = largest_divisor(numPPRanksRequested);
781 /* Check if the largest divisor is more than numPPRanks ^ (2/3) */
782 if (largestDivisor * largestDivisor * largestDivisor > numPPRanksRequested * numPPRanksRequested)
785 "The number of ranks selected for particle-particle work (%d) "
786 "contains a large prime factor %d. In most cases this will lead to "
787 "bad performance. Choose a number with smaller prime factors or "
788 "set the decomposition (option -dd) manually.",
789 numPPRanksRequested, largestDivisor);
794 /*! \brief Return the number of PME-only ranks used by the simulation
796 * If the user did not choose a number, then decide for them. */
797 static int getNumPmeOnlyRanksToUse(const gmx::MDLogger& mdlog,
798 const DomdecOptions& options,
799 const gmx_mtop_t& mtop,
800 const t_inputrec& ir,
802 const int numRanksRequested)
805 const char* extraMessage = "";
807 if (options.numCells[XX] > 0)
809 if (options.numPmeRanks >= 0)
811 // TODO mdrun should give a fatal error with a non-PME input file and -npme > 0
812 numPmeOnlyRanks = options.numPmeRanks;
816 // When the DD grid is set explicitly and -npme is set to auto,
817 // don't use PME ranks. We check later if the DD grid is
818 // compatible with the total number of ranks.
824 if (!EEL_PME(ir.coulombtype))
826 // TODO mdrun should give a fatal error with a non-PME input file and -npme > 0
831 if (options.numPmeRanks >= 0)
833 numPmeOnlyRanks = options.numPmeRanks;
837 // Controls the automated choice of when to use separate PME-only ranks.
838 const int minRankCountToDefaultToSeparatePmeRanks = 19;
840 if (numRanksRequested < minRankCountToDefaultToSeparatePmeRanks)
844 ", as there are too few total\n"
845 " ranks for efficient splitting";
849 numPmeOnlyRanks = guess_npme(mdlog, mtop, ir, box, numRanksRequested);
850 extraMessage = ", as guessed by mdrun";
855 GMX_RELEASE_ASSERT(numPmeOnlyRanks <= numRanksRequested,
856 "Cannot have more PME ranks than total ranks");
857 if (EEL_PME(ir.coulombtype))
859 GMX_LOG(mdlog.info).appendTextFormatted("Using %d separate PME ranks%s", numPmeOnlyRanks, extraMessage);
862 return numPmeOnlyRanks;
865 /*! \brief Sets the order of the DD dimensions, returns the number of DD dimensions */
866 static int set_dd_dim(const gmx::IVec& numDDCells, const DDSettings& ddSettings, ivec* dims)
869 if (ddSettings.useDDOrderZYX)
871 /* Decomposition order z,y,x */
872 for (int dim = DIM - 1; dim >= 0; dim--)
874 if (numDDCells[dim] > 1)
876 (*dims)[ndim++] = dim;
882 /* Decomposition order x,y,z */
883 for (int dim = 0; dim < DIM; dim++)
885 if (numDDCells[dim] > 1)
887 (*dims)[ndim++] = dim;
894 /* Set dim[0] to avoid extra checks on ndim in several places */
901 DDGridSetup getDDGridSetup(const gmx::MDLogger& mdlog,
903 const int numRanksRequested,
904 const DomdecOptions& options,
905 const DDSettings& ddSettings,
906 const DDSystemInfo& systemInfo,
907 const real cellSizeLimit,
908 const gmx_mtop_t& mtop,
909 const t_inputrec& ir,
911 gmx::ArrayRef<const gmx::RVec> xGlobal,
914 int numPmeOnlyRanks = getNumPmeOnlyRanksToUse(mdlog, options, mtop, ir, box, numRanksRequested);
916 if (ddSettings.request1DAnd1Pulse && (numRanksRequested - numPmeOnlyRanks == 1))
918 // With only one PP rank, there will not be a need for
919 // GPU-based halo exchange that wants to request that any DD
920 // has only 1 dimension and 1 pulse.
921 return DDGridSetup{};
924 gmx::IVec numDomains;
925 if (options.numCells[XX] > 0)
927 numDomains = gmx::IVec(options.numCells);
928 const ivec numDomainsLegacyIvec = { numDomains[XX], numDomains[YY], numDomains[ZZ] };
929 set_ddbox_cr(*cr, &numDomainsLegacyIvec, ir, box, xGlobal, ddbox);
933 set_ddbox_cr(*cr, nullptr, ir, box, xGlobal, ddbox);
937 numDomains = optimizeDDCells(mdlog, numRanksRequested, numPmeOnlyRanks, cellSizeLimit,
938 ddSettings.request1DAnd1Pulse, mtop, box, *ddbox, ir, systemInfo);
942 /* Communicate the information set by the master to all ranks */
943 gmx_bcast(sizeof(numDomains), numDomains, cr);
944 if (EEL_PME(ir.coulombtype))
946 gmx_bcast(sizeof(numPmeOnlyRanks), &numPmeOnlyRanks, cr);
949 DDGridSetup ddGridSetup;
950 ddGridSetup.numPmeOnlyRanks = numPmeOnlyRanks;
951 ddGridSetup.numDomains[XX] = numDomains[XX];
952 ddGridSetup.numDomains[YY] = numDomains[YY];
953 ddGridSetup.numDomains[ZZ] = numDomains[ZZ];
954 ddGridSetup.numDDDimensions = set_dd_dim(numDomains, ddSettings, &ddGridSetup.ddDimensions);