2 * This file is part of the GROMACS molecular simulation package.
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.
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.
43 #include "gromacs/legacyheaders/domdec.h"
44 #include "gromacs/legacyheaders/names.h"
45 #include "gromacs/legacyheaders/network.h"
46 #include "gromacs/legacyheaders/perf_est.h"
47 #include "gromacs/legacyheaders/typedefs.h"
48 #include "gromacs/legacyheaders/types/commrec.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/utility/smalloc.h"
52 /* Margin for setting up the DD grid */
53 #define DD_GRID_MARGIN_PRES_SCALE 1.05
55 static int factorize(int n, int **fac, int **mfac)
61 gmx_fatal(FARGS, "Can only factorize positive integers.");
64 /* Decompose n in factors */
73 if (ndiv == 0 || (*fac)[ndiv-1] != d)
87 static gmx_bool largest_divisor(int n)
89 int ndiv, *div, *mdiv, ldiv;
91 ndiv = factorize(n, &div, &mdiv);
99 static int lcd(int n1, int n2)
104 for (i = 2; (i <= n1 && i <= n2); i++)
106 if (n1 % i == 0 && n2 % i == 0)
115 static gmx_bool fits_pme_ratio(int nnodes, int npme, float ratio)
117 return ((double)npme/(double)nnodes > 0.95*ratio);
120 static gmx_bool fits_pp_pme_perf(int nnodes, int npme, float ratio)
122 int ndiv, *div, *mdiv, ldiv;
123 int npp_root3, npme_root2;
125 ndiv = factorize(nnodes-npme, &div, &mdiv);
130 npp_root3 = static_cast<int>(std::pow(nnodes-npme, 1.0/3.0) + 0.5);
131 npme_root2 = static_cast<int>(std::sqrt(static_cast<double>(npme)) + 0.5);
133 /* The check below gives a reasonable division:
134 * factor 5 allowed at 5 or more PP nodes,
135 * factor 7 allowed at 49 or more PP nodes.
137 if (ldiv > 3 + npp_root3)
142 /* Check if the number of PP and PME nodes have a reasonable sized
143 * denominator in common, such that we can use 2D PME decomposition
144 * when required (which requires nx_pp == nx_pme).
145 * The factor of 2 allows for a maximum ratio of 2^2=4
146 * between nx_pme and ny_pme.
148 if (lcd(nnodes-npme, npme)*2 < npme_root2)
153 /* Does this division gives a reasonable PME load? */
154 return fits_pme_ratio(nnodes, npme, ratio);
157 static int guess_npme(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir, matrix box,
163 ratio = pme_load_estimate(mtop, ir, box);
167 fprintf(fplog, "Guess for relative PME load: %.2f\n", ratio);
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.
175 if (!fits_pme_ratio(nnodes, nnodes/2, ratio))
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.
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.
189 npme = (nnodes + 15)/16;
190 while (npme <= nnodes/3)
192 if (nnodes % npme == 0)
194 /* Note that fits_perf might change the PME grid,
195 * in the current implementation it does not.
197 if (fits_pp_pme_perf(nnodes, npme, ratio))
206 /* Try any possible number for npme */
208 while (npme <= nnodes/2)
210 /* Note that fits_perf may change the PME grid */
211 if (fits_pp_pme_perf(nnodes, npme, ratio))
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 */
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",
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",
244 static int div_up(int n, int f)
246 return (n + f - 1)/f;
249 real comm_box_frac(ivec dd_nc, real cutoff, gmx_ddbox_t *ddbox)
255 for (i = 0; i < DIM; i++)
257 bt[i] = ddbox->box_size[i]*ddbox->skew_fac[i];
258 nw[i] = dd_nc[i]*cutoff/bt[i];
262 for (i = 0; i < DIM; i++)
267 for (j = i+1; j < DIM; j++)
271 comm_vol += nw[i]*nw[j]*M_PI/4;
272 for (k = j+1; k < DIM; k++)
276 comm_vol += nw[i]*nw[j]*nw[k]*M_PI/6;
287 static gmx_bool inhomogeneous_z(const t_inputrec *ir)
289 return ((EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD) &&
290 ir->ePBC == epbcXYZ && ir->ewald_geometry == eewg3DC);
293 /* Avoid integer overflows */
294 static float comm_pme_cost_vol(int npme, int a, int b, int c)
300 comm_vol *= div_up(a, npme);
301 comm_vol *= div_up(b, npme);
306 static float comm_cost_est(real limit, real cutoff,
307 matrix box, gmx_ddbox_t *ddbox,
308 int natoms, t_inputrec *ir,
310 int npme_tot, ivec nc)
312 ivec npme = {1, 1, 1};
313 int i, j, nk, overlap;
315 float comm_vol, comm_vol_xf, comm_pme, cost_pbcdx;
316 /* This is the cost of a pbc_dx call relative to the cost
317 * of communicating the coordinate and force of an atom.
318 * This will be machine dependent.
319 * These factors are for x86 with SMP or Infiniband.
321 float pbcdx_rect_fac = 0.1;
322 float pbcdx_tric_fac = 0.2;
325 /* Check the DD algorithm restrictions */
326 if ((ir->ePBC == epbcXY && ir->nwall < 2 && nc[ZZ] > 1) ||
327 (ir->ePBC == epbcSCREW && (nc[XX] == 1 || nc[YY] > 1 || nc[ZZ] > 1)))
332 if (inhomogeneous_z(ir) && nc[ZZ] > 1)
337 assert(ddbox->npbcdim <= DIM);
339 /* Check if the triclinic requirements are met */
340 for (i = 0; i < DIM; i++)
342 for (j = i+1; j < ddbox->npbcdim; j++)
344 if (box[j][i] != 0 || ir->deform[j][i] != 0 ||
345 (ir->epc != epcNO && ir->compress[j][i] != 0))
347 if (nc[j] > 1 && nc[i] == 1)
355 for (i = 0; i < DIM; i++)
357 bt[i] = ddbox->box_size[i]*ddbox->skew_fac[i];
359 /* Without PBC and with 2 cells, there are no lower limits on the cell size */
360 if (!(i >= ddbox->npbcdim && nc[i] <= 2) && bt[i] < nc[i]*limit)
364 /* With PBC, check if the cut-off fits in nc[i]-1 cells */
365 if (i < ddbox->npbcdim && nc[i] > 1 && (nc[i] - 1)*bt[i] < nc[i]*cutoff)
373 /* The following choices should match those
374 * in init_domain_decomposition in domdec.c.
376 if (nc[XX] == 1 && nc[YY] > 1)
381 else if (nc[YY] == 1)
388 /* Will we use 1D or 2D PME decomposition? */
389 npme[XX] = (npme_tot % nc[XX] == 0) ? nc[XX] : npme_tot;
390 npme[YY] = npme_tot/npme[XX];
394 /* When two dimensions are (nearly) equal, use more cells
395 * for the smallest index, so the decomposition does not
396 * depend sensitively on the rounding of the box elements.
398 for (i = 0; i < DIM; i++)
400 for (j = i+1; j < DIM; j++)
402 /* Check if the box size is nearly identical,
403 * in that case we prefer nx > ny and ny > nz.
405 if (fabs(bt[j] - bt[i]) < 0.01*bt[i] && nc[j] > nc[i])
407 /* The XX/YY check is a bit compact. If nc[YY]==npme[YY]
408 * this means the swapped nc has nc[XX]==npme[XX],
409 * and we can also swap X and Y for PME.
411 /* Check if dimension i and j are equivalent for PME.
412 * For x/y: if nc[YY]!=npme[YY], we can not swap x/y
413 * For y/z: we can not have PME decomposition in z
416 !((i == XX && j == YY && nc[YY] != npme[YY]) ||
417 (i == YY && j == ZZ && npme[YY] > 1)))
425 /* This function determines only half of the communication cost.
426 * All PP, PME and PP-PME communication is symmetric
427 * and the "back"-communication cost is identical to the forward cost.
430 comm_vol = comm_box_frac(nc, cutoff, ddbox);
433 for (i = 0; i < 2; i++)
435 /* Determine the largest volume for PME x/f redistribution */
436 if (nc[i] % npme[i] != 0)
440 comm_vol_xf = (npme[i] == 2 ? 1.0/3.0 : 0.5);
444 comm_vol_xf = 1.0 - lcd(nc[i], npme[i])/(double)npme[i];
446 comm_pme += 3*natoms*comm_vol_xf;
449 /* Grid overlap communication */
452 nk = (i == 0 ? ir->nkx : ir->nky);
453 overlap = (nk % npme[i] == 0 ? ir->pme_order-1 : ir->pme_order);
461 /* Old line comm_pme += npme[i]*overlap*ir->nkx*ir->nky*ir->nkz/nk; */
465 /* PME FFT communication volume.
466 * This only takes the communication into account and not imbalance
467 * in the calculation. But the imbalance in communication and calculation
468 * are similar and therefore these formulas also prefer load balance
469 * in the FFT and pme_solve calculation.
471 comm_pme += comm_pme_cost_vol(npme[YY], ir->nky, ir->nkz, ir->nkx);
472 comm_pme += comm_pme_cost_vol(npme[XX], ir->nkx, ir->nky, ir->nkz);
474 /* Add cost of pbc_dx for bondeds */
476 if ((nc[XX] == 1 || nc[YY] == 1) || (nc[ZZ] == 1 && ir->ePBC != epbcXY))
478 if ((ddbox->tric_dir[XX] && nc[XX] == 1) ||
479 (ddbox->tric_dir[YY] && nc[YY] == 1))
481 cost_pbcdx = pbcdxr*pbcdx_tric_fac;
485 cost_pbcdx = pbcdxr*pbcdx_rect_fac;
492 "nc %2d %2d %2d %2d %2d vol pp %6.4f pbcdx %6.4f pme %9.3e tot %9.3e\n",
493 nc[XX], nc[YY], nc[ZZ], npme[XX], npme[YY],
494 comm_vol, cost_pbcdx, comm_pme,
495 3*natoms*(comm_vol + cost_pbcdx) + comm_pme);
498 return 3*natoms*(comm_vol + cost_pbcdx) + comm_pme;
501 static void assign_factors(gmx_domdec_t *dd,
502 real limit, real cutoff,
503 matrix box, gmx_ddbox_t *ddbox,
504 int natoms, t_inputrec *ir,
505 float pbcdxr, int npme,
506 int ndiv, int *div, int *mdiv, ivec ir_try, ivec opt)
513 ce = comm_cost_est(limit, cutoff, box, ddbox,
514 natoms, ir, pbcdxr, npme, ir_try);
515 if (ce >= 0 && (opt[XX] == 0 ||
516 ce < comm_cost_est(limit, cutoff, box, ddbox,
520 copy_ivec(ir_try, opt);
526 for (x = mdiv[0]; x >= 0; x--)
528 for (i = 0; i < x; i++)
530 ir_try[XX] *= div[0];
532 for (y = mdiv[0]-x; y >= 0; y--)
534 for (i = 0; i < y; i++)
536 ir_try[YY] *= div[0];
538 for (i = 0; i < mdiv[0]-x-y; i++)
540 ir_try[ZZ] *= div[0];
544 assign_factors(dd, limit, cutoff, box, ddbox, natoms, ir, pbcdxr, npme,
545 ndiv-1, div+1, mdiv+1, ir_try, opt);
547 for (i = 0; i < mdiv[0]-x-y; i++)
549 ir_try[ZZ] /= div[0];
551 for (i = 0; i < y; i++)
553 ir_try[YY] /= div[0];
556 for (i = 0; i < x; i++)
558 ir_try[XX] /= div[0];
563 static real optimize_ncells(FILE *fplog,
564 int nnodes_tot, int npme_only,
565 gmx_bool bDynLoadBal, real dlb_scale,
566 gmx_mtop_t *mtop, matrix box, gmx_ddbox_t *ddbox,
569 real cellsize_limit, real cutoff,
570 gmx_bool bInterCGBondeds,
573 int npp, npme, ndiv, *div, *mdiv, d, nmax;
574 gmx_bool bExcl_pbcdx;
579 limit = cellsize_limit;
585 npp = nnodes_tot - npme_only;
586 if (EEL_PME(ir->coulombtype))
588 npme = (npme_only > 0 ? npme_only : npp);
597 /* For Ewald exclusions pbc_dx is not called */
599 (IR_EXCL_FORCES(*ir) && !EEL_FULL(ir->coulombtype));
600 pbcdxr = (double)n_bonded_dx(mtop, bExcl_pbcdx)/(double)mtop->natoms;
604 /* Every molecule is a single charge group: no pbc required */
607 /* Add a margin for DLB and/or pressure scaling */
610 if (dlb_scale >= 1.0)
612 gmx_fatal(FARGS, "The value for option -dds should be smaller than 1");
616 fprintf(fplog, "Scaling the initial minimum size with 1/%g (option -dds) = %g\n", dlb_scale, 1/dlb_scale);
620 else if (ir->epc != epcNO)
624 fprintf(fplog, "To account for pressure scaling, scaling the initial minimum size with %g\n", DD_GRID_MARGIN_PRES_SCALE);
625 limit *= DD_GRID_MARGIN_PRES_SCALE;
631 fprintf(fplog, "Optimizing the DD grid for %d cells with a minimum initial size of %.3f nm\n", npp, limit);
633 if (inhomogeneous_z(ir))
635 fprintf(fplog, "Ewald_geometry=%s: assuming inhomogeneous particle distribution in z, will not decompose in z.\n", eewg_names[ir->ewald_geometry]);
640 fprintf(fplog, "The maximum allowed number of cells is:");
641 for (d = 0; d < DIM; d++)
643 nmax = (int)(ddbox->box_size[d]*ddbox->skew_fac[d]/limit);
644 if (d >= ddbox->npbcdim && nmax < 2)
648 if (d == ZZ && inhomogeneous_z(ir))
652 fprintf(fplog, " %c %d", 'X' + d, nmax);
654 fprintf(fplog, "\n");
660 fprintf(debug, "Average nr of pbc_dx calls per atom %.2f\n", pbcdxr);
663 /* Decompose npp in factors */
664 ndiv = factorize(npp, &div, &mdiv);
670 assign_factors(dd, limit, cutoff, box, ddbox, mtop->natoms, ir, pbcdxr,
671 npme, ndiv, div, mdiv, itry, nc);
679 real dd_choose_grid(FILE *fplog,
680 t_commrec *cr, gmx_domdec_t *dd, t_inputrec *ir,
681 gmx_mtop_t *mtop, matrix box, gmx_ddbox_t *ddbox,
682 gmx_bool bDynLoadBal, real dlb_scale,
683 real cellsize_limit, real cutoff_dd,
684 gmx_bool bInterCGBondeds)
686 gmx_int64_t nnodes_div, ldiv;
691 nnodes_div = cr->nnodes;
692 if (EEL_PME(ir->coulombtype))
694 if (cr->npmenodes > 0)
699 "Cannot have separate PME ranks with 2 or fewer ranks");
701 if (cr->npmenodes >= cr->nnodes)
704 "Cannot have %d separate PME ranks with just %d total ranks",
705 cr->npmenodes, cr->nnodes);
708 /* If the user purposely selected the number of PME nodes,
709 * only check for large primes in the PP node count.
711 nnodes_div -= cr->npmenodes;
721 ldiv = largest_divisor(nnodes_div);
722 /* Check if the largest divisor is more than nnodes^2/3 */
723 if (ldiv*ldiv*ldiv > nnodes_div*nnodes_div)
725 gmx_fatal(FARGS, "The number of ranks you selected (%d) contains a large prime factor %d. In most cases this will lead to bad performance. Choose a number with smaller prime factors or set the decomposition (option -dd) manually.",
730 if (EEL_PME(ir->coulombtype))
732 if (cr->npmenodes < 0)
734 /* Use PME nodes when the number of nodes is more than 16 */
735 if (cr->nnodes <= 18)
740 fprintf(fplog, "Using %d separate PME ranks, as there are too few total\n ranks for efficient splitting\n", cr->npmenodes);
745 cr->npmenodes = guess_npme(fplog, mtop, ir, box, cr->nnodes);
748 fprintf(fplog, "Using %d separate PME ranks, as guessed by mdrun\n", cr->npmenodes);
756 fprintf(fplog, "Using %d separate PME ranks, per user request\n", cr->npmenodes);
761 limit = optimize_ncells(fplog, cr->nnodes, cr->npmenodes,
762 bDynLoadBal, dlb_scale,
763 mtop, box, ddbox, ir, dd,
764 cellsize_limit, cutoff_dd,
772 /* Communicate the information set by the master to all nodes */
773 gmx_bcast(sizeof(dd->nc), dd->nc, cr);
774 if (EEL_PME(ir->coulombtype))
776 gmx_bcast(sizeof(ir->nkx), &ir->nkx, cr);
777 gmx_bcast(sizeof(ir->nky), &ir->nky, cr);
778 gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);