2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2008, The GROMACS development team.
6 * Copyright (c) 2012,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gromacs/math/functions.h"
44 #include "gromacs/math/utilities.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/mdtypes/commrec.h"
47 #include "gromacs/mdtypes/inputrec.h"
48 #include "gromacs/mdtypes/md_enums.h"
49 #include "gromacs/nbnxm/nbnxm.h"
50 #include "gromacs/simd/simd.h"
51 #include "gromacs/topology/ifunc.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/fatalerror.h"
55 /* Computational cost of bonded, non-bonded and PME calculations.
56 * This will be machine dependent.
57 * The numbers are only used for estimating the relative cost of PME vs PP,
58 * so only relative numbers matter.
59 * The numbers here are accurate cycle counts for Haswell in single precision
60 * compiled with gcc5.2. A correction factor for other architectures is given
61 * by simd_cycle_factor().
62 * In double precision PME mesh is slightly cheaper, although not so much
63 * that the numbers need to be adjusted.
66 /* Cost of a pair interaction in the "group" cut-off scheme */
67 static const double c_group_fq = 18.0;
68 static const double c_group_qlj_cut = 18.0;
69 static const double c_group_qlj_tab = 24.0;
70 static const double c_group_lj_cut = 12.0;
71 static const double c_group_lj_tab = 21.0;
72 /* Cost of 1 water with one Q/LJ atom */
73 static const double c_group_qljw_cut = 24.0;
74 static const double c_group_qljw_tab = 27.0;
75 /* Cost of 1 water with one Q atom or with 1/3 water (LJ negligible) */
76 static const double c_group_qw = 21.0;
78 /* Cost of a pair interaction in the "Verlet" cut-off scheme, QEXP is Ewald */
79 static const double c_nbnxn_lj = 2.5;
80 static const double c_nbnxn_qrf_lj = 2.9;
81 static const double c_nbnxn_qrf = 2.4;
82 static const double c_nbnxn_qexp_lj = 4.2;
83 static const double c_nbnxn_qexp = 3.8;
84 /* Extra cost for expensive LJ interaction, e.g. pot-switch or LJ-PME */
85 static const double c_nbnxn_ljexp_add = 1.0;
87 /* Cost of the different components of PME. */
88 /* Cost of particle reordering and redistribution (no SIMD correction).
89 * This will be zero without MPI and can be very high with load imbalance.
90 * Thus we use an approximate value for medium parallelization.
92 static const double c_pme_redist = 100.0;
93 /* Cost of q spreading and force interpolation per charge. This part almost
94 * doesn't accelerate with SIMD, so we don't use SIMD correction.
96 static const double c_pme_spread = 5.0;
97 /* Cost of fft's, will be multiplied with 2 N log2(N) (no SIMD correction)
98 * Without MPI the number is 2-3, depending on grid factors and thread count.
99 * We take the high limit to be on the safe side and account for some MPI
100 * communication cost, which will dominate at high parallelization.
102 static const double c_pme_fft = 3.0;
103 /* Cost of pme_solve, will be multiplied with N */
104 static const double c_pme_solve = 9.0;
106 /* Cost of a bonded interaction divided by the number of distances calculations
107 * required in one interaction. The actual cost is nearly propotional to this.
109 static const double c_bond = 25.0;
112 #if GMX_SIMD_HAVE_REAL
113 static const gmx_bool bHaveSIMD = TRUE;
115 static const gmx_bool bHaveSIMD = FALSE;
118 /* Gives a correction factor for the currently compiled SIMD implementations
119 * versus the reference used for the coefficients above (8-wide SIMD with FMA).
120 * bUseSIMD sets if we asking for plain-C (FALSE) or SIMD (TRUE) code.
122 static double simd_cycle_factor(gmx_bool bUseSIMD)
124 /* The (average) ratio of the time taken by plain-C force calculations
125 * relative to SIMD versions, for the reference platform Haswell:
126 * 8-wide SIMD with FMA, factor: sqrt(2*8)*1.25 = 5.
127 * This factor is used for normalization in simd_cycle_factor().
129 const double simd_cycle_no_simd = 5.0;
132 #if GMX_SIMD_HAVE_REAL
135 /* We never get full speed-up of a factor GMX_SIMD_REAL_WIDTH.
136 * The actual speed-up depends very much on gather+scatter overhead,
137 * which is different for different bonded and non-bonded kernels.
138 * As a rough, but actually not bad, approximation we use a sqrt
139 * dependence on the width which gives a factor 4 for width=8.
141 speedup = std::sqrt(2.0*GMX_SIMD_REAL_WIDTH);
142 #if GMX_SIMD_HAVE_FMA
143 /* FMA tends to give a bit more speedup */
154 gmx_incons("gmx_cycle_factor() compiled without SIMD called with bUseSIMD=TRUE");
156 /* No SIMD, no speedup */
160 /* Return speed compared to the reference (Haswell).
161 * For x86 SIMD, the nbnxn kernels are relatively much slower on
162 * Sandy/Ivy Bridge than Haswell, but that only leads to a too high
163 * PME load estimate on SB/IB, which is erring on the safe side.
165 return simd_cycle_no_simd/speedup;
168 void count_bonded_distances(const gmx_mtop_t *mtop, const t_inputrec *ir,
169 double *ndistance_c, double *ndistance_simd)
172 double nonsimd_step_frac;
174 double ndtot_c, ndtot_simd;
175 #if GMX_SIMD_HAVE_REAL
176 gmx_bool bSimdBondeds = TRUE;
178 gmx_bool bSimdBondeds = FALSE;
181 bExcl = (ir->cutoff_scheme == ecutsGROUP && inputrecExclForces(ir)
182 && !EEL_FULL(ir->coulombtype));
186 /* We only have SIMD versions of these bondeds without energy and
187 * without shift-forces, we take that into account here.
189 if (ir->nstcalcenergy > 0)
191 nonsimd_step_frac = 1.0/ir->nstcalcenergy;
195 nonsimd_step_frac = 0;
197 if (ir->epc != epcNO && 1.0/ir->nstpcouple > nonsimd_step_frac)
199 nonsimd_step_frac = 1.0/ir->nstpcouple;
204 nonsimd_step_frac = 1;
207 /* Count the number of pbc_rvec_sub calls required for bonded interactions.
208 * This number is also roughly proportional to the computational cost.
212 for (const gmx_molblock_t &molb : mtop->molblock)
214 const gmx_moltype_t *molt = &mtop->moltype[molb.type];
215 for (ftype = 0; ftype < F_NRE; ftype++)
219 if (interaction_function[ftype].flags & IF_BOND)
221 double nd_c, nd_simd;
225 /* For all interactions, except for the three exceptions
226 * in the switch below, #distances = #atoms - 1.
236 /* These bonded potentially use SIMD */
241 nd_c = nonsimd_step_frac *(NRAL(ftype) - 1);
242 nd_simd = (1 - nonsimd_step_frac)*(NRAL(ftype) - 1);
245 nd_c = NRAL(ftype) - 1;
248 nbonds = molb.nmol*molt->ilist[ftype].size()/(1 + NRAL(ftype));
249 ndtot_c += nbonds*nd_c;
250 ndtot_simd += nbonds*nd_simd;
255 ndtot_c += molb.nmol*(molt->excls.nra - molt->atoms.nr)/2.;
261 fprintf(debug, "nr. of distance calculations in bondeds: C %.1f SIMD %.1f\n", ndtot_c, ndtot_simd);
264 if (ndistance_c != nullptr)
266 *ndistance_c = ndtot_c;
268 if (ndistance_simd != nullptr)
270 *ndistance_simd = ndtot_simd;
274 static void pp_group_load(const gmx_mtop_t *mtop, const t_inputrec *ir,
276 int *nq_tot, int *nlj_tot,
278 gmx_bool *bChargePerturbed, gmx_bool *bTypePerturbed)
280 int atnr, cg, a0, ncqlj, ncq, nclj;
281 gmx_bool bBHAM, bLJcut, bWater, bQ, bLJ;
282 int nw, nqlj, nq, nlj;
283 double fq, fqlj, flj, fqljw, fqw;
285 bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
287 bLJcut = ((ir->vdwtype == evdwCUT) && !bBHAM);
289 /* Computational cost of bonded, non-bonded and PME calculations.
290 * This will be machine dependent.
291 * The numbers here are accurate for Intel Core2 and AMD Athlon 64
292 * in single precision. In double precision PME mesh is slightly cheaper,
293 * although not so much that the numbers need to be adjusted.
296 fqlj = (bLJcut ? c_group_qlj_cut : c_group_qlj_tab);
297 flj = (bLJcut ? c_group_lj_cut : c_group_lj_tab);
298 /* Cost of 1 water with one Q/LJ atom */
299 fqljw = (bLJcut ? c_group_qljw_cut : c_group_qljw_tab);
300 /* Cost of 1 water with one Q atom or with 1/3 water (LJ negligible) */
303 gmx::ArrayRef<const t_iparams> iparams = mtop->ffparams.iparams;
304 atnr = mtop->ffparams.atnr;
309 *bChargePerturbed = FALSE;
310 *bTypePerturbed = FALSE;
311 for (const gmx_molblock_t &molb : mtop->molblock)
313 const gmx_moltype_t *molt = &mtop->moltype[molb.type];
314 const t_atom *atom = molt->atoms.atom;
316 for (cg = 0; cg < molt->cgs.nr; cg++)
323 while (a < molt->cgs.index[cg+1])
325 bQ = (atom[a].q != 0 || atom[a].qB != 0);
326 bLJ = (iparams[(atnr+1)*atom[a].type].lj.c6 != 0 ||
327 iparams[(atnr+1)*atom[a].type].lj.c12 != 0);
328 if (atom[a].q != atom[a].qB)
330 *bChargePerturbed = TRUE;
332 if (atom[a].type != atom[a].typeB)
334 *bTypePerturbed = TRUE;
336 /* This if this atom fits into water optimization */
337 if (!((a == a0 && bQ && bLJ) ||
338 (a == a0+1 && bQ && !bLJ) ||
339 (a == a0+2 && bQ && !bLJ && atom[a].q == atom[a-1].q) ||
340 (a == a0+3 && !bQ && bLJ)))
367 nqlj += molb.nmol*ncqlj;
369 nlj += molb.nmol*nclj;
374 *nq_tot = nq + nqlj + nw*3;
375 *nlj_tot = nlj + nqlj + nw;
379 fprintf(debug, "nw %d nqlj %d nq %d nlj %d\n", nw, nqlj, nq, nlj);
382 /* For the PP non-bonded cost it is (unrealistically) assumed
383 * that all atoms are distributed homogeneously in space.
384 * Factor 3 is used because a water molecule has 3 atoms
385 * (and TIP4P effectively has 3 interactions with (water) atoms)).
387 *cost_pp = 0.5*(fqljw*nw*nqlj +
388 fqw *nw*(3*nw + nq) +
390 fq *nq*(3*nw + nqlj + nq) +
391 flj *nlj*(nw + nqlj + nlj))
392 *4/3*M_PI*ir->rlist*ir->rlist*ir->rlist/det(box);
394 *cost_pp *= simd_cycle_factor(bHaveSIMD);
397 static void pp_verlet_load(const gmx_mtop_t *mtop, const t_inputrec *ir,
399 int *nq_tot, int *nlj_tot,
401 gmx_bool *bChargePerturbed, gmx_bool *bTypePerturbed)
403 int atnr, a, nqlj, nq, nlj;
406 double c_qlj, c_q, c_lj;
409 /* Conversion factor for reference vs SIMD kernel performance.
410 * The factor is about right for SSE2/4, but should be 2 higher for AVX256.
413 const real nbnxn_refkernel_fac = 4.0;
415 const real nbnxn_refkernel_fac = 8.0;
418 bQRF = (EEL_RF(ir->coulombtype) || ir->coulombtype == eelCUT);
420 gmx::ArrayRef<const t_iparams> iparams = mtop->ffparams.iparams;
421 atnr = mtop->ffparams.atnr;
424 *bChargePerturbed = FALSE;
425 *bTypePerturbed = FALSE;
426 for (const gmx_molblock_t &molb : mtop->molblock)
428 const gmx_moltype_t *molt = &mtop->moltype[molb.type];
429 const t_atom *atom = molt->atoms.atom;
430 for (a = 0; a < molt->atoms.nr; a++)
432 if (atom[a].q != 0 || atom[a].qB != 0)
434 if (iparams[(atnr+1)*atom[a].type].lj.c6 != 0 ||
435 iparams[(atnr+1)*atom[a].type].lj.c12 != 0)
444 if (atom[a].q != atom[a].qB)
446 *bChargePerturbed = TRUE;
448 if (atom[a].type != atom[a].typeB)
450 *bTypePerturbed = TRUE;
455 nlj = mtop->natoms - nqlj - nq;
458 *nlj_tot = nqlj + nlj;
460 /* Effective cut-off for cluster pair list of 4x4 or 4x8 atoms.
461 * This choice should match the one of pick_nbnxn_kernel_cpu().
462 * TODO: Make this function use pick_nbnxn_kernel_cpu().
464 #if GMX_SIMD_HAVE_REAL && ((GMX_SIMD_REAL_WIDTH == 8 && defined GMX_SIMD_HAVE_FMA) || GMX_SIMD_REAL_WIDTH > 8)
469 r_eff = ir->rlist + nbnxn_get_rlist_effective_inc(j_cluster_size, mtop->natoms/det(box));
471 /* The average number of pairs per atom */
472 nppa = 0.5*4/3*M_PI*r_eff*r_eff*r_eff*mtop->natoms/det(box);
476 fprintf(debug, "nqlj %d nq %d nlj %d rlist %.3f r_eff %.3f pairs per atom %.1f\n",
477 nqlj, nq, nlj, ir->rlist, r_eff, nppa);
480 /* Determine the cost per pair interaction */
481 c_qlj = (bQRF ? c_nbnxn_qrf_lj : c_nbnxn_qexp_lj);
482 c_q = (bQRF ? c_nbnxn_qrf : c_nbnxn_qexp);
484 if (ir->vdw_modifier == eintmodPOTSWITCH || EVDW_PME(ir->vdwtype))
486 c_qlj += c_nbnxn_ljexp_add;
487 c_lj += c_nbnxn_ljexp_add;
489 if (EVDW_PME(ir->vdwtype) && ir->ljpme_combination_rule == eljpmeLB)
491 /* We don't have LJ-PME LB comb. rule kernels, we use slow kernels */
492 c_qlj *= nbnxn_refkernel_fac;
493 c_q *= nbnxn_refkernel_fac;
494 c_lj *= nbnxn_refkernel_fac;
497 /* For the PP non-bonded cost it is (unrealistically) assumed
498 * that all atoms are distributed homogeneously in space.
500 *cost_pp = (nqlj*c_qlj + nq*c_q + nlj*c_lj)*nppa;
502 *cost_pp *= simd_cycle_factor(bHaveSIMD);
505 float pme_load_estimate(const gmx_mtop_t *mtop, const t_inputrec *ir,
509 gmx_bool bChargePerturbed, bTypePerturbed;
510 double ndistance_c, ndistance_simd;
511 double cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve, cost_pme;
514 /* Computational cost of bonded, non-bonded and PME calculations.
515 * This will be machine dependent.
516 * The numbers here are accurate for Intel Core2 and AMD Athlon 64
517 * in single precision. In double precision PME mesh is slightly cheaper,
518 * although not so much that the numbers need to be adjusted.
521 count_bonded_distances(mtop, ir, &ndistance_c, &ndistance_simd);
522 /* C_BOND is the cost for bonded interactions with SIMD implementations,
523 * so we need to scale the number of bonded interactions for which there
524 * are only C implementations to the number of SIMD equivalents.
526 cost_bond = c_bond*(ndistance_c *simd_cycle_factor(FALSE) +
527 ndistance_simd*simd_cycle_factor(bHaveSIMD));
529 if (ir->cutoff_scheme == ecutsGROUP)
531 pp_group_load(mtop, ir, box,
532 &nq_tot, &nlj_tot, &cost_pp,
533 &bChargePerturbed, &bTypePerturbed);
537 pp_verlet_load(mtop, ir, box,
538 &nq_tot, &nlj_tot, &cost_pp,
539 &bChargePerturbed, &bTypePerturbed);
547 int gridNkzFactor = int{
550 if (EEL_PME(ir->coulombtype))
552 double grid = ir->nkx*ir->nky*gridNkzFactor;
554 int f = ((ir->efep != efepNO && bChargePerturbed) ? 2 : 1);
555 cost_redist += c_pme_redist*nq_tot;
556 cost_spread += f*c_pme_spread*nq_tot*gmx::power3(ir->pme_order);
557 cost_fft += f*c_pme_fft*grid*std::log(grid)/std::log(2.0);
558 cost_solve += f*c_pme_solve*grid*simd_cycle_factor(bHaveSIMD);
561 if (EVDW_PME(ir->vdwtype))
563 double grid = ir->nkx*ir->nky*gridNkzFactor;
565 int f = ((ir->efep != efepNO && bTypePerturbed) ? 2 : 1);
566 if (ir->ljpme_combination_rule == eljpmeLB)
568 /* LB combination rule: we have 7 mesh terms */
571 cost_redist += c_pme_redist*nlj_tot;
572 cost_spread += f*c_pme_spread*nlj_tot*gmx::power3(ir->pme_order);
573 cost_fft += f*c_pme_fft*2*grid*std::log(grid)/std::log(2.0);
574 cost_solve += f*c_pme_solve*grid*simd_cycle_factor(bHaveSIMD);
577 cost_pme = cost_redist + cost_spread + cost_fft + cost_solve;
579 ratio = cost_pme/(cost_bond + cost_pp + cost_pme);
590 cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve);
592 fprintf(debug, "Estimate for relative PME load: %.3f\n", ratio);