1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This file is part of Gromacs Copyright (c) 1991-2008
5 * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
7 * This program is free software; you can redistribute it and/or
8 * modify it under the terms of the GNU General Public License
9 * as published by the Free Software Foundation; either version 2
10 * of the License, or (at your option) any later version.
12 * To help us fund GROMACS development, we humbly ask that you cite
13 * the research papers on the package. Check out http://www.gromacs.org
16 * Gnomes, ROck Monsters And Chili Sauce
34 /* Margin for setting up the DD grid */
35 #define DD_GRID_MARGIN_PRES_SCALE 1.05
37 static int factorize(int n,int **fac,int **mfac)
43 gmx_fatal(FARGS, "Can only factorize positive integers.");
46 /* Decompose n in factors */
55 if (ndiv == 0 || (*fac)[ndiv-1] != d)
69 static gmx_bool largest_divisor(int n)
71 int ndiv,*div,*mdiv,ldiv;
73 ndiv = factorize(n,&div,&mdiv);
81 static int lcd(int n1,int n2)
86 for(i=2; (i<=n1 && i<=n2); i++)
88 if (n1 % i == 0 && n2 % i == 0)
97 static gmx_bool fits_pme_ratio(int nnodes,int npme,float ratio)
99 return ((double)npme/(double)nnodes > 0.95*ratio);
102 static gmx_bool fits_pp_pme_perf(FILE *fplog,
103 t_inputrec *ir,matrix box,gmx_mtop_t *mtop,
104 int nnodes,int npme,float ratio)
106 int ndiv,*div,*mdiv,ldiv;
107 int npp_root3,npme_root2;
109 ndiv = factorize(nnodes-npme,&div,&mdiv);
114 npp_root3 = (int)(pow(nnodes-npme,1.0/3.0) + 0.5);
115 npme_root2 = (int)(sqrt(npme) + 0.5);
117 /* The check below gives a reasonable division:
118 * factor 5 allowed at 5 or more PP nodes,
119 * factor 7 allowed at 49 or more PP nodes.
121 if (ldiv > 3 + npp_root3)
126 /* Check if the number of PP and PME nodes have a reasonable sized
127 * denominator in common, such that we can use 2D PME decomposition
128 * when required (which requires nx_pp == nx_pme).
129 * The factor of 2 allows for a maximum ratio of 2^2=4
130 * between nx_pme and ny_pme.
132 if (lcd(nnodes-npme,npme)*2 < npme_root2)
137 /* Does this division gives a reasonable PME load? */
138 return fits_pme_ratio(nnodes,npme,ratio);
141 static int guess_npme(FILE *fplog,gmx_mtop_t *mtop,t_inputrec *ir,matrix box,
148 ratio = pme_load_estimate(mtop,ir,box);
152 fprintf(fplog,"Guess for relative PME load: %.2f\n",ratio);
155 /* We assume the optimal node ratio is close to the load ratio.
156 * The communication load is neglected,
157 * but (hopefully) this will balance out between PP and PME.
160 if (!fits_pme_ratio(nnodes,nnodes/2,ratio))
162 /* We would need more than nnodes/2 PME only nodes,
163 * which is not possible. Since the PME load is very high,
164 * we will not loose much performance when all nodes do PME.
170 /* First try to find npme as a factor of nnodes up to nnodes/3.
171 * We start with a minimum PME node fraction of 1/16
172 * and avoid ratios which lead to large prime factors in nnodes-npme.
174 npme = (nnodes + 15)/16;
175 while (npme <= nnodes/3) {
176 if (nnodes % npme == 0)
178 /* Note that fits_perf might change the PME grid,
179 * in the current implementation it does not.
181 if (fits_pp_pme_perf(fplog,ir,box,mtop,nnodes,npme,ratio))
190 /* Try any possible number for npme */
192 while (npme <= nnodes/2)
194 /* Note that fits_perf may change the PME grid */
195 if (fits_pp_pme_perf(fplog,ir,box,mtop,nnodes,npme,ratio))
204 gmx_fatal(FARGS,"Could not find an appropriate number of separate PME nodes. i.e. >= %5f*#nodes (%d) and <= #nodes/2 (%d) and reasonable performance wise (grid_x=%d, grid_y=%d).\n"
205 "Use the -npme option of mdrun or change the number of processors or the PME grid dimensions, see the manual for details.",
206 ratio,(int)(0.95*ratio*nnodes+0.5),nnodes/2,ir->nkx,ir->nky);
207 /* Keep the compiler happy */
215 "Will use %d particle-particle and %d PME only nodes\n"
216 "This is a guess, check the performance at the end of the log file\n",
220 "Will use %d particle-particle and %d PME only nodes\n"
221 "This is a guess, check the performance at the end of the log file\n",
228 static int div_up(int n,int f)
230 return (n + f - 1)/f;
233 real comm_box_frac(ivec dd_nc,real cutoff,gmx_ddbox_t *ddbox)
241 bt[i] = ddbox->box_size[i]*ddbox->skew_fac[i];
242 nw[i] = dd_nc[i]*cutoff/bt[i];
253 for(j=i+1; j<DIM; j++)
257 comm_vol += nw[i]*nw[j]*M_PI/4;
258 for(k=j+1; k<DIM; k++)
262 comm_vol += nw[i]*nw[j]*nw[k]*M_PI/6;
273 static gmx_bool inhomogeneous_z(const t_inputrec *ir)
275 return ((EEL_PME(ir->coulombtype) || ir->coulombtype==eelEWALD) &&
276 ir->ePBC==epbcXYZ && ir->ewald_geometry==eewg3DC);
279 /* Avoid integer overflows */
280 static float comm_pme_cost_vol(int npme, int a, int b, int c)
286 comm_vol *= div_up(a, npme);
287 comm_vol *= div_up(b, npme);
292 static float comm_cost_est(gmx_domdec_t *dd,real limit,real cutoff,
293 matrix box,gmx_ddbox_t *ddbox,
294 int natoms,t_inputrec *ir,
296 int npme_tot,ivec nc)
299 int i,j,k,nk,overlap;
301 float comm_vol,comm_vol_xf,comm_pme,cost_pbcdx;
302 /* This is the cost of a pbc_dx call relative to the cost
303 * of communicating the coordinate and force of an atom.
304 * This will be machine dependent.
305 * These factors are for x86 with SMP or Infiniband.
307 float pbcdx_rect_fac = 0.1;
308 float pbcdx_tric_fac = 0.2;
311 /* Check the DD algorithm restrictions */
312 if ((ir->ePBC == epbcXY && ir->nwall < 2 && nc[ZZ] > 1) ||
313 (ir->ePBC == epbcSCREW && (nc[XX] == 1 || nc[YY] > 1 || nc[ZZ] > 1)))
318 if (inhomogeneous_z(ir) && nc[ZZ] > 1)
323 assert(ddbox->npbcdim<=DIM);
325 /* Check if the triclinic requirements are met */
328 for(j=i+1; j<ddbox->npbcdim; j++)
330 if (box[j][i] != 0 || ir->deform[j][i] != 0 ||
331 (ir->epc != epcNO && ir->compress[j][i] != 0))
333 if (nc[j] > 1 && nc[i] == 1)
343 bt[i] = ddbox->box_size[i]*ddbox->skew_fac[i];
345 /* Without PBC there are no cell size limits with 2 cells */
346 if (!(i >= ddbox->npbcdim && nc[i] <= 2) && bt[i] < nc[i]*limit)
354 /* The following choices should match those
355 * in init_domain_decomposition in domdec.c.
357 if (nc[XX] == 1 && nc[YY] > 1)
362 else if (nc[YY] == 1)
369 /* Will we use 1D or 2D PME decomposition? */
370 npme[XX] = (npme_tot % nc[XX] == 0) ? nc[XX] : npme_tot;
371 npme[YY] = npme_tot/npme[XX];
375 /* When two dimensions are (nearly) equal, use more cells
376 * for the smallest index, so the decomposition does not
377 * depend sensitively on the rounding of the box elements.
381 for(j=i+1; j<DIM; j++)
383 /* Check if the box size is nearly identical,
384 * in that case we prefer nx > ny and ny > nz.
386 if (fabs(bt[j] - bt[i]) < 0.01*bt[i] && nc[j] > nc[i])
388 /* The XX/YY check is a bit compact. If nc[YY]==npme[YY]
389 * this means the swapped nc has nc[XX]==npme[XX],
390 * and we can also swap X and Y for PME.
392 /* Check if dimension i and j are equivalent for PME.
393 * For x/y: if nc[YY]!=npme[YY], we can not swap x/y
394 * For y/z: we can not have PME decomposition in z
397 !((i == XX && j == YY && nc[YY] != npme[YY]) ||
398 (i == YY && j == ZZ && npme[YY] > 1)))
406 /* This function determines only half of the communication cost.
407 * All PP, PME and PP-PME communication is symmetric
408 * and the "back"-communication cost is identical to the forward cost.
411 comm_vol = comm_box_frac(nc,cutoff,ddbox);
416 /* Determine the largest volume for PME x/f redistribution */
417 if (nc[i] % npme[i] != 0)
421 comm_vol_xf = (npme[i]==2 ? 1.0/3.0 : 0.5);
425 comm_vol_xf = 1.0 - lcd(nc[i],npme[i])/(double)npme[i];
427 comm_pme += 3*natoms*comm_vol_xf;
430 /* Grid overlap communication */
433 nk = (i==0 ? ir->nkx : ir->nky);
434 overlap = (nk % npme[i] == 0 ? ir->pme_order-1 : ir->pme_order);
442 /* Old line comm_pme += npme[i]*overlap*ir->nkx*ir->nky*ir->nkz/nk; */
446 /* PME FFT communication volume.
447 * This only takes the communication into account and not imbalance
448 * in the calculation. But the imbalance in communication and calculation
449 * are similar and therefore these formulas also prefer load balance
450 * in the FFT and pme_solve calculation.
452 comm_pme += comm_pme_cost_vol(npme[YY], ir->nky, ir->nkz, ir->nkx);
453 comm_pme += comm_pme_cost_vol(npme[XX], ir->nkx, ir->nky, ir->nkz);
455 /* Add cost of pbc_dx for bondeds */
457 if ((nc[XX] == 1 || nc[YY] == 1) || (nc[ZZ] == 1 && ir->ePBC != epbcXY))
459 if ((ddbox->tric_dir[XX] && nc[XX] == 1) ||
460 (ddbox->tric_dir[YY] && nc[YY] == 1))
462 cost_pbcdx = pbcdxr*pbcdx_tric_fac;
466 cost_pbcdx = pbcdxr*pbcdx_rect_fac;
473 "nc %2d %2d %2d %2d %2d vol pp %6.4f pbcdx %6.4f pme %9.3e tot %9.3e\n",
474 nc[XX],nc[YY],nc[ZZ],npme[XX],npme[YY],
475 comm_vol,cost_pbcdx,comm_pme,
476 3*natoms*(comm_vol + cost_pbcdx) + comm_pme);
479 return 3*natoms*(comm_vol + cost_pbcdx) + comm_pme;
482 static void assign_factors(gmx_domdec_t *dd,
483 real limit,real cutoff,
484 matrix box,gmx_ddbox_t *ddbox,
485 int natoms,t_inputrec *ir,
486 float pbcdxr,int npme,
487 int ndiv,int *div,int *mdiv,ivec ir_try,ivec opt)
494 ce = comm_cost_est(dd,limit,cutoff,box,ddbox,
495 natoms,ir,pbcdxr,npme,ir_try);
496 if (ce >= 0 && (opt[XX] == 0 ||
497 ce < comm_cost_est(dd,limit,cutoff,box,ddbox,
501 copy_ivec(ir_try,opt);
507 for(x=mdiv[0]; x>=0; x--)
511 ir_try[XX] *= div[0];
513 for(y=mdiv[0]-x; y>=0; y--)
517 ir_try[YY] *= div[0];
519 for(i=0; i<mdiv[0]-x-y; i++)
521 ir_try[ZZ] *= div[0];
525 assign_factors(dd,limit,cutoff,box,ddbox,natoms,ir,pbcdxr,npme,
526 ndiv-1,div+1,mdiv+1,ir_try,opt);
528 for(i=0; i<mdiv[0]-x-y; i++)
530 ir_try[ZZ] /= div[0];
534 ir_try[YY] /= div[0];
539 ir_try[XX] /= div[0];
544 static real optimize_ncells(FILE *fplog,
545 int nnodes_tot,int npme_only,
546 gmx_bool bDynLoadBal,real dlb_scale,
547 gmx_mtop_t *mtop,matrix box,gmx_ddbox_t *ddbox,
550 real cellsize_limit,real cutoff,
551 gmx_bool bInterCGBondeds,gmx_bool bInterCGMultiBody,
554 int npp,npme,ndiv,*div,*mdiv,d,nmax;
555 gmx_bool bExcl_pbcdx;
560 limit = cellsize_limit;
566 npp = nnodes_tot - npme_only;
567 if (EEL_PME(ir->coulombtype))
569 npme = (npme_only > 0 ? npme_only : npp);
578 /* For Ewald exclusions pbc_dx is not called */
580 (IR_EXCL_FORCES(*ir) && !EEL_FULL(ir->coulombtype));
581 pbcdxr = (double)n_bonded_dx(mtop,bExcl_pbcdx)/(double)mtop->natoms;
585 /* Every molecule is a single charge group: no pbc required */
588 /* Add a margin for DLB and/or pressure scaling */
591 if (dlb_scale >= 1.0)
593 gmx_fatal(FARGS,"The value for option -dds should be smaller than 1");
597 fprintf(fplog,"Scaling the initial minimum size with 1/%g (option -dds) = %g\n",dlb_scale,1/dlb_scale);
601 else if (ir->epc != epcNO)
605 fprintf(fplog,"To account for pressure scaling, scaling the initial minimum size with %g\n",DD_GRID_MARGIN_PRES_SCALE);
606 limit *= DD_GRID_MARGIN_PRES_SCALE;
612 fprintf(fplog,"Optimizing the DD grid for %d cells with a minimum initial size of %.3f nm\n",npp,limit);
614 if (inhomogeneous_z(ir))
616 fprintf(fplog,"Ewald_geometry=%s: assuming inhomogeneous particle distribution in z, will not decompose in z.\n",eewg_names[ir->ewald_geometry]);
621 fprintf(fplog,"The maximum allowed number of cells is:");
624 nmax = (int)(ddbox->box_size[d]*ddbox->skew_fac[d]/limit);
625 if (d >= ddbox->npbcdim && nmax < 2)
629 if (d == ZZ && inhomogeneous_z(ir))
633 fprintf(fplog," %c %d",'X' + d,nmax);
641 fprintf(debug,"Average nr of pbc_dx calls per atom %.2f\n",pbcdxr);
644 /* Decompose npp in factors */
645 ndiv = factorize(npp,&div,&mdiv);
651 assign_factors(dd,limit,cutoff,box,ddbox,mtop->natoms,ir,pbcdxr,
652 npme,ndiv,div,mdiv,itry,nc);
660 real dd_choose_grid(FILE *fplog,
661 t_commrec *cr,gmx_domdec_t *dd,t_inputrec *ir,
662 gmx_mtop_t *mtop,matrix box,gmx_ddbox_t *ddbox,
663 gmx_bool bDynLoadBal,real dlb_scale,
664 real cellsize_limit,real cutoff_dd,
665 gmx_bool bInterCGBondeds,gmx_bool bInterCGMultiBody)
667 gmx_large_int_t nnodes_div,ldiv;
672 nnodes_div = cr->nnodes;
673 if (EEL_PME(ir->coulombtype))
675 if (cr->npmenodes > 0)
680 "Can not have separate PME nodes with 2 or less nodes");
682 if (cr->npmenodes >= cr->nnodes)
685 "Can not have %d separate PME nodes with just %d total nodes",
686 cr->npmenodes, cr->nnodes);
689 /* If the user purposely selected the number of PME nodes,
690 * only check for large primes in the PP node count.
692 nnodes_div -= cr->npmenodes;
702 ldiv = largest_divisor(nnodes_div);
703 /* Check if the largest divisor is more than nnodes^2/3 */
704 if (ldiv*ldiv*ldiv > nnodes_div*nnodes_div)
706 gmx_fatal(FARGS,"The number of nodes 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.",
711 if (EEL_PME(ir->coulombtype))
713 if (cr->npmenodes < 0)
715 /* Use PME nodes when the number of nodes is more than 16 */
716 if (cr->nnodes <= 18)
721 fprintf(fplog,"Using %d separate PME nodes, as there are too few total\n nodes for efficient splitting\n",cr->npmenodes);
726 cr->npmenodes = guess_npme(fplog,mtop,ir,box,cr->nnodes);
729 fprintf(fplog,"Using %d separate PME nodes, as guessed by mdrun\n",cr->npmenodes);
737 fprintf(fplog,"Using %d separate PME nodes, per user request\n",cr->npmenodes);
742 limit = optimize_ncells(fplog,cr->nnodes,cr->npmenodes,
743 bDynLoadBal,dlb_scale,
744 mtop,box,ddbox,ir,dd,
745 cellsize_limit,cutoff_dd,
746 bInterCGBondeds,bInterCGMultiBody,
753 /* Communicate the information set by the master to all nodes */
754 gmx_bcast(sizeof(dd->nc),dd->nc,cr);
755 if (EEL_PME(ir->coulombtype))
757 gmx_bcast(sizeof(ir->nkx),&ir->nkx,cr);
758 gmx_bcast(sizeof(ir->nky),&ir->nky,cr);
759 gmx_bcast(sizeof(cr->npmenodes),&cr->npmenodes,cr);