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, 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/vec.h"
45 #include "gromacs/mdlib/nbnxn_consts.h"
46 #include "gromacs/mdlib/nbnxn_search.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/md_enums.h"
49 #include "gromacs/simd/simd.h"
50 #include "gromacs/topology/ifunc.h"
51 #include "gromacs/topology/topology.h"
52 #include "gromacs/utility/fatalerror.h"
54 /* Computational cost of bonded, non-bonded and PME calculations.
55 * This will be machine dependent.
56 * The numbers are only used for estimating the relative cost of PME vs PP,
57 * so only relative numbers matter.
58 * The numbers here are accurate cycle counts for Haswell in single precision
59 * compiled with gcc5.2. A correction factor for other architectures is given
60 * by simd_cycle_factor().
61 * In double precision PME mesh is slightly cheaper, although not so much
62 * that the numbers need to be adjusted.
65 /* Cost of a pair interaction in the "group" cut-off scheme */
66 static const double c_group_fq = 18.0;
67 static const double c_group_qlj_cut = 18.0;
68 static const double c_group_qlj_tab = 24.0;
69 static const double c_group_lj_cut = 12.0;
70 static const double c_group_lj_tab = 21.0;
71 /* Cost of 1 water with one Q/LJ atom */
72 static const double c_group_qljw_cut = 24.0;
73 static const double c_group_qljw_tab = 27.0;
74 /* Cost of 1 water with one Q atom or with 1/3 water (LJ negligible) */
75 static const double c_group_qw = 21.0;
77 /* Cost of a pair interaction in the "Verlet" cut-off scheme, QEXP is Ewald */
78 static const double c_nbnxn_lj = 2.5;
79 static const double c_nbnxn_qrf_lj = 2.9;
80 static const double c_nbnxn_qrf = 2.4;
81 static const double c_nbnxn_qexp_lj = 4.2;
82 static const double c_nbnxn_qexp = 3.8;
83 /* Extra cost for expensive LJ interaction, e.g. pot-switch or LJ-PME */
84 static const double c_nbnxn_ljexp_add = 1.0;
86 /* Cost of the different components of PME. */
87 /* Cost of particle reordering and redistribution (no SIMD correction).
88 * This will be zero without MPI and can be very high with load imbalance.
89 * Thus we use an approximate value for medium parallelization.
91 static const double c_pme_redist = 100.0;
92 /* Cost of q spreading and force interpolation per charge. This part almost
93 * doesn't accelerate with SIMD, so we don't use SIMD correction.
95 static const double c_pme_spread = 5.0;
96 /* Cost of fft's, will be multiplied with 2 N log2(N) (no SIMD correction)
97 * Without MPI the number is 2-3, depending on grid factors and thread count.
98 * We take the high limit to be on the safe side and account for some MPI
99 * communication cost, which will dominate at high parallelization.
101 static const double c_pme_fft = 3.0;
102 /* Cost of pme_solve, will be multiplied with N */
103 static const double c_pme_solve = 9.0;
105 /* Cost of a bonded interaction divided by the number of distances calculations
106 * required in one interaction. The actual cost is nearly propotional to this.
108 static const double c_bond = 25.0;
111 #if GMX_SIMD_HAVE_REAL
112 static const gmx_bool bHaveSIMD = TRUE;
114 static const gmx_bool bHaveSIMD = FALSE;
117 /* Gives a correction factor for the currently compiled SIMD implementations
118 * versus the reference used for the coefficients above (8-wide SIMD with FMA).
119 * bUseSIMD sets if we asking for plain-C (FALSE) or SIMD (TRUE) code.
121 static double simd_cycle_factor(gmx_bool bUseSIMD)
123 /* The (average) ratio of the time taken by plain-C force calculations
124 * relative to SIMD versions, for the reference platform Haswell:
125 * 8-wide SIMD with FMA, factor: sqrt(2*8)*1.25 = 5.
126 * This factor is used for normalization in simd_cycle_factor().
128 const double simd_cycle_no_simd = 5.0;
131 #if GMX_SIMD_HAVE_REAL
134 /* We never get full speed-up of a factor GMX_SIMD_REAL_WIDTH.
135 * The actual speed-up depends very much on gather+scatter overhead,
136 * which is different for different bonded and non-bonded kernels.
137 * As a rough, but actually not bad, approximation we use a sqrt
138 * dependence on the width which gives a factor 4 for width=8.
140 speedup = std::sqrt(2.0*GMX_SIMD_REAL_WIDTH);
141 #if GMX_SIMD_HAVE_FMA
142 /* FMA tends to give a bit more speedup */
153 gmx_incons("gmx_cycle_factor() compiled without SIMD called with bUseSIMD=TRUE");
155 /* No SIMD, no speedup */
159 /* Return speed compared to the reference (Haswell).
160 * For x86 SIMD, the nbnxn kernels are relatively much slower on
161 * Sandy/Ivy Bridge than Haswell, but that only leads to a too high
162 * PME load estimate on SB/IB, which is erring on the safe side.
164 return simd_cycle_no_simd/speedup;
167 void count_bonded_distances(const gmx_mtop_t *mtop, const t_inputrec *ir,
168 double *ndistance_c, double *ndistance_simd)
171 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 #if defined _ICC && __ICC == 1400 || defined __ICL && __ICL == 1400
213 #pragma novector /* Work-around for incorrect vectorization */
215 for (mb = 0; mb < mtop->nmolblock; mb++)
217 molt = &mtop->moltype[mtop->molblock[mb].type];
218 nmol = mtop->molblock[mb].nmol;
219 for (ftype = 0; ftype < F_NRE; ftype++)
223 if (interaction_function[ftype].flags & IF_BOND)
225 double nd_c, nd_simd;
229 /* For all interactions, except for the three exceptions
230 * in the switch below, #distances = #atoms - 1.
240 /* These bonded potentially use SIMD */
244 nd_c = nonsimd_step_frac *(NRAL(ftype) - 1);
245 nd_simd = (1 - nonsimd_step_frac)*(NRAL(ftype) - 1);
248 nd_c = NRAL(ftype) - 1;
251 nbonds = nmol*molt->ilist[ftype].nr/(1 + NRAL(ftype));
252 ndtot_c += nbonds*nd_c;
253 ndtot_simd += nbonds*nd_simd;
258 ndtot_c += nmol*(molt->excls.nra - molt->atoms.nr)/2;
264 fprintf(debug, "nr. of distance calculations in bondeds: C %.1f SIMD %.1f\n", ndtot_c, ndtot_simd);
267 if (ndistance_c != NULL)
269 *ndistance_c = ndtot_c;
271 if (ndistance_simd != NULL)
273 *ndistance_simd = ndtot_simd;
277 static void pp_group_load(const gmx_mtop_t *mtop, const t_inputrec *ir,
279 int *nq_tot, int *nlj_tot,
281 gmx_bool *bChargePerturbed, gmx_bool *bTypePerturbed)
284 int mb, nmol, atnr, cg, a, a0, ncqlj, ncq, nclj;
285 gmx_bool bBHAM, bLJcut, bWater, bQ, bLJ;
286 int nw, nqlj, nq, nlj;
287 double fq, fqlj, flj, fqljw, fqw;
291 bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
293 bLJcut = ((ir->vdwtype == evdwCUT) && !bBHAM);
295 /* Computational cost of bonded, non-bonded and PME calculations.
296 * This will be machine dependent.
297 * The numbers here are accurate for Intel Core2 and AMD Athlon 64
298 * in single precision. In double precision PME mesh is slightly cheaper,
299 * although not so much that the numbers need to be adjusted.
302 fqlj = (bLJcut ? c_group_qlj_cut : c_group_qlj_tab);
303 flj = (bLJcut ? c_group_lj_cut : c_group_lj_tab);
304 /* Cost of 1 water with one Q/LJ atom */
305 fqljw = (bLJcut ? c_group_qljw_cut : c_group_qljw_tab);
306 /* Cost of 1 water with one Q atom or with 1/3 water (LJ negligible) */
309 iparams = mtop->ffparams.iparams;
310 atnr = mtop->ffparams.atnr;
315 *bChargePerturbed = FALSE;
316 *bTypePerturbed = FALSE;
317 for (mb = 0; mb < mtop->nmolblock; mb++)
319 molt = &mtop->moltype[mtop->molblock[mb].type];
320 atom = molt->atoms.atom;
321 nmol = mtop->molblock[mb].nmol;
323 for (cg = 0; cg < molt->cgs.nr; cg++)
330 while (a < molt->cgs.index[cg+1])
332 bQ = (atom[a].q != 0 || atom[a].qB != 0);
333 bLJ = (iparams[(atnr+1)*atom[a].type].lj.c6 != 0 ||
334 iparams[(atnr+1)*atom[a].type].lj.c12 != 0);
335 if (atom[a].q != atom[a].qB)
337 *bChargePerturbed = TRUE;
339 if (atom[a].type != atom[a].typeB)
341 *bTypePerturbed = TRUE;
343 /* This if this atom fits into water optimization */
344 if (!((a == a0 && bQ && bLJ) ||
345 (a == a0+1 && bQ && !bLJ) ||
346 (a == a0+2 && bQ && !bLJ && atom[a].q == atom[a-1].q) ||
347 (a == a0+3 && !bQ && bLJ)))
381 *nq_tot = nq + nqlj + nw*3;
382 *nlj_tot = nlj + nqlj + nw;
386 fprintf(debug, "nw %d nqlj %d nq %d nlj %d\n", nw, nqlj, nq, nlj);
389 /* For the PP non-bonded cost it is (unrealistically) assumed
390 * that all atoms are distributed homogeneously in space.
391 * Factor 3 is used because a water molecule has 3 atoms
392 * (and TIP4P effectively has 3 interactions with (water) atoms)).
394 *cost_pp = 0.5*(fqljw*nw*nqlj +
395 fqw *nw*(3*nw + nq) +
397 fq *nq*(3*nw + nqlj + nq) +
398 flj *nlj*(nw + nqlj + nlj))
399 *4/3*M_PI*ir->rlist*ir->rlist*ir->rlist/det(box);
401 *cost_pp *= simd_cycle_factor(bHaveSIMD);
404 static void pp_verlet_load(const gmx_mtop_t *mtop, const t_inputrec *ir,
406 int *nq_tot, int *nlj_tot,
408 gmx_bool *bChargePerturbed, gmx_bool *bTypePerturbed)
411 int mb, nmol, atnr, a, nqlj, nq, nlj;
416 double c_qlj, c_q, c_lj;
419 /* Conversion factor for reference vs SIMD kernel performance.
420 * The factor is about right for SSE2/4, but should be 2 higher for AVX256.
423 const real nbnxn_refkernel_fac = 4.0;
425 const real nbnxn_refkernel_fac = 8.0;
428 bQRF = (EEL_RF(ir->coulombtype) || ir->coulombtype == eelCUT);
430 iparams = mtop->ffparams.iparams;
431 atnr = mtop->ffparams.atnr;
434 *bChargePerturbed = FALSE;
435 *bTypePerturbed = FALSE;
436 for (mb = 0; mb < mtop->nmolblock; mb++)
438 molt = &mtop->moltype[mtop->molblock[mb].type];
439 atom = molt->atoms.atom;
440 nmol = mtop->molblock[mb].nmol;
441 for (a = 0; a < molt->atoms.nr; a++)
443 if (atom[a].q != 0 || atom[a].qB != 0)
445 if (iparams[(atnr+1)*atom[a].type].lj.c6 != 0 ||
446 iparams[(atnr+1)*atom[a].type].lj.c12 != 0)
455 if (atom[a].q != atom[a].qB)
457 *bChargePerturbed = TRUE;
459 if (atom[a].type != atom[a].typeB)
461 *bTypePerturbed = TRUE;
466 nlj = mtop->natoms - nqlj - nq;
469 *nlj_tot = nqlj + nlj;
471 /* Effective cut-off for cluster pair list of 4x4 or 4x8 atoms.
472 * This choice should match the one of pick_nbnxn_kernel_cpu().
473 * TODO: Make this function use pick_nbnxn_kernel_cpu().
475 #if GMX_SIMD_HAVE_REAL && ((GMX_SIMD_REAL_WIDTH == 8 && defined GMX_SIMD_HAVE_FMA) || GMX_SIMD_REAL_WIDTH > 8)
480 r_eff = ir->rlist + nbnxn_get_rlist_effective_inc(j_cluster_size, mtop->natoms/det(box));
482 /* The average number of pairs per atom */
483 nppa = 0.5*4/3*M_PI*r_eff*r_eff*r_eff*mtop->natoms/det(box);
487 fprintf(debug, "nqlj %d nq %d nlj %d rlist %.3f r_eff %.3f pairs per atom %.1f\n",
488 nqlj, nq, nlj, ir->rlist, r_eff, nppa);
491 /* Determine the cost per pair interaction */
492 c_qlj = (bQRF ? c_nbnxn_qrf_lj : c_nbnxn_qexp_lj);
493 c_q = (bQRF ? c_nbnxn_qrf : c_nbnxn_qexp);
495 if (ir->vdw_modifier == eintmodPOTSWITCH || EVDW_PME(ir->vdwtype))
497 c_qlj += c_nbnxn_ljexp_add;
498 c_lj += c_nbnxn_ljexp_add;
500 if (EVDW_PME(ir->vdwtype) && ir->ljpme_combination_rule == eljpmeLB)
502 /* We don't have LJ-PME LB comb. rule kernels, we use slow kernels */
503 c_qlj *= nbnxn_refkernel_fac;
504 c_q *= nbnxn_refkernel_fac;
505 c_lj *= nbnxn_refkernel_fac;
508 /* For the PP non-bonded cost it is (unrealistically) assumed
509 * that all atoms are distributed homogeneously in space.
511 *cost_pp = (nqlj*c_qlj + nq*c_q + nlj*c_lj)*nppa;
513 *cost_pp *= simd_cycle_factor(bHaveSIMD);
516 float pme_load_estimate(const gmx_mtop_t *mtop, const t_inputrec *ir,
520 gmx_bool bChargePerturbed, bTypePerturbed;
521 double ndistance_c, ndistance_simd;
522 double cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve, cost_pme;
525 /* Computational cost of bonded, non-bonded and PME calculations.
526 * This will be machine dependent.
527 * The numbers here are accurate for Intel Core2 and AMD Athlon 64
528 * in single precision. In double precision PME mesh is slightly cheaper,
529 * although not so much that the numbers need to be adjusted.
532 count_bonded_distances(mtop, ir, &ndistance_c, &ndistance_simd);
533 /* C_BOND is the cost for bonded interactions with SIMD implementations,
534 * so we need to scale the number of bonded interactions for which there
535 * are only C implementations to the number of SIMD equivalents.
537 cost_bond = c_bond*(ndistance_c *simd_cycle_factor(FALSE) +
538 ndistance_simd*simd_cycle_factor(bHaveSIMD));
540 if (ir->cutoff_scheme == ecutsGROUP)
542 pp_group_load(mtop, ir, box,
543 &nq_tot, &nlj_tot, &cost_pp,
544 &bChargePerturbed, &bTypePerturbed);
548 pp_verlet_load(mtop, ir, box,
549 &nq_tot, &nlj_tot, &cost_pp,
550 &bChargePerturbed, &bTypePerturbed);
558 if (EEL_PME(ir->coulombtype))
560 double grid = ir->nkx*ir->nky*((ir->nkz + 1)/2);
562 int f = ((ir->efep != efepNO && bChargePerturbed) ? 2 : 1);
563 cost_redist += c_pme_redist*nq_tot;
564 cost_spread += f*c_pme_spread*nq_tot*gmx::power3(ir->pme_order);
565 cost_fft += f*c_pme_fft*grid*std::log(grid)/std::log(2.0);
566 cost_solve += f*c_pme_solve*grid*simd_cycle_factor(bHaveSIMD);
569 if (EVDW_PME(ir->vdwtype))
571 double grid = ir->nkx*ir->nky*((ir->nkz + 1)/2);
573 int f = ((ir->efep != efepNO && bTypePerturbed) ? 2 : 1);
574 if (ir->ljpme_combination_rule == eljpmeLB)
576 /* LB combination rule: we have 7 mesh terms */
579 cost_redist += c_pme_redist*nlj_tot;
580 cost_spread += f*c_pme_spread*nlj_tot*gmx::power3(ir->pme_order);
581 cost_fft += f*c_pme_fft*2*grid*std::log(grid)/std::log(2.0);
582 cost_solve += f*c_pme_solve*grid*simd_cycle_factor(bHaveSIMD);
585 cost_pme = cost_redist + cost_spread + cost_fft + cost_solve;
587 ratio = cost_pme/(cost_bond + cost_pp + cost_pme);
598 cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve);
600 fprintf(debug, "Estimate for relative PME load: %.3f\n", ratio);