2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
36 * \brief Common functions for the different NBNXN GPU implementations.
38 * \author Berk Hess <hess@kth.se>
40 * \ingroup module_nbnxm
45 #include "gromacs/domdec/domdec.h"
46 #include "gromacs/domdec/domdec_struct.h"
47 #include "gromacs/hardware/hw_info.h"
48 #include "gromacs/mdlib/gmx_omp_nthreads.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/mdtypes/forcerec.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/nbnxm/atomdata.h"
53 #include "gromacs/nbnxm/gpu_data_mgmt.h"
54 #include "gromacs/nbnxm/nbnxm.h"
55 #include "gromacs/nbnxm/nbnxm_geometry.h"
56 #include "gromacs/nbnxm/nbnxm_simd.h"
57 #include "gromacs/nbnxm/pairlist_tuning.h"
58 #include "gromacs/nbnxm/pairlistset.h"
59 #include "gromacs/simd/simd.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/logger.h"
69 /*! \brief Returns whether CPU SIMD support exists for the given inputrec
71 * If the return value is FALSE and fplog/cr != NULL, prints a fallback
72 * message to fplog/stderr.
74 static gmx_bool nbnxn_simd_supported(const gmx::MDLogger &mdlog,
77 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
79 /* LJ PME with LB combination rule does 7 mesh operations.
80 * This so slow that we don't compile SIMD non-bonded kernels
82 GMX_LOG(mdlog.warning).asParagraph().appendText("LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels");
89 /*! \brief Returns the most suitable CPU kernel type and Ewald handling */
90 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
93 const gmx_hw_info_t gmx_unused &hardwareInfo)
95 *kernel_type = nbnxnk4x4_PlainC;
96 *ewald_excl = ewaldexclTable;
100 #ifdef GMX_NBNXN_SIMD_4XN
101 *kernel_type = nbnxnk4xN_SIMD_4xN;
103 #ifdef GMX_NBNXN_SIMD_2XNN
104 *kernel_type = nbnxnk4xN_SIMD_2xNN;
107 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
108 /* We need to choose if we want 2x(N+N) or 4xN kernels.
109 * This is based on the SIMD acceleration choice and CPU information
110 * detected at runtime.
112 * 4xN calculates more (zero) interactions, but has less pair-search
113 * work and much better kernel instruction scheduling.
115 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
116 * which doesn't have FMA, both the analytical and tabulated Ewald
117 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
118 * 2x(4+4) because it results in significantly fewer pairs.
119 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
120 * 10% with HT, 50% without HT. As we currently don't detect the actual
121 * use of HT, use 4x8 to avoid a potential performance hit.
122 * On Intel Haswell 4x8 is always faster.
124 *kernel_type = nbnxnk4xN_SIMD_4xN;
126 #if !GMX_SIMD_HAVE_FMA
127 if (EEL_PME_EWALD(ir->coulombtype) ||
128 EVDW_PME(ir->vdwtype))
130 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
131 * There are enough instructions to make 2x(4+4) efficient.
133 *kernel_type = nbnxnk4xN_SIMD_2xNN;
136 if (hardwareInfo.haveAmdZenCpu)
138 /* One 256-bit FMA per cycle makes 2xNN faster */
139 *kernel_type = nbnxnk4xN_SIMD_2xNN;
141 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
144 if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
146 #ifdef GMX_NBNXN_SIMD_4XN
147 *kernel_type = nbnxnk4xN_SIMD_4xN;
149 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
152 if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
154 #ifdef GMX_NBNXN_SIMD_2XNN
155 *kernel_type = nbnxnk4xN_SIMD_2xNN;
157 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
161 /* Analytical Ewald exclusion correction is only an option in
163 * Since table lookup's don't parallelize with SIMD, analytical
164 * will probably always be faster for a SIMD width of 8 or more.
165 * With FMA analytical is sometimes faster for a width if 4 as well.
166 * In single precision, this is faster on Bulldozer.
168 #if GMX_SIMD_REAL_WIDTH >= 8 || \
169 (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)
170 /* On AMD Zen, tabulated Ewald kernels are faster on all 4 combinations
171 * of single or double precision and 128 or 256-bit AVX2.
173 if (!hardwareInfo.haveAmdZenCpu)
175 *ewald_excl = ewaldexclAnalytical;
178 if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
180 *ewald_excl = ewaldexclTable;
182 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
184 *ewald_excl = ewaldexclAnalytical;
191 const char *lookup_kernel_name(int kernel_type)
193 const char *returnvalue = nullptr;
197 returnvalue = "not set";
199 case nbnxnk4x4_PlainC:
200 returnvalue = "plain C";
202 case nbnxnk4xN_SIMD_4xN:
203 case nbnxnk4xN_SIMD_2xNN:
205 returnvalue = "SIMD";
207 returnvalue = "not available";
210 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
211 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
215 gmx_fatal(FARGS, "Illegal kernel type selected");
220 /*! \brief Returns the most suitable kernel type and Ewald handling */
221 static void pick_nbnxn_kernel(const gmx::MDLogger &mdlog,
222 gmx_bool use_simd_kernels,
223 const gmx_hw_info_t &hardwareInfo,
225 EmulateGpuNonbonded emulateGpu,
226 const t_inputrec *ir,
229 gmx_bool bDoNonbonded)
231 GMX_RELEASE_ASSERT(kernel_type, "Need a valid kernel_type pointer");
233 *kernel_type = nbnxnkNotSet;
234 *ewald_excl = ewaldexclTable;
236 if (emulateGpu == EmulateGpuNonbonded::Yes)
238 *kernel_type = nbnxnk8x8x8_PlainC;
242 GMX_LOG(mdlog.warning).asParagraph().appendText("Emulating a GPU run on the CPU (slow)");
247 *kernel_type = nbnxnk8x8x8_GPU;
250 if (*kernel_type == nbnxnkNotSet)
252 if (use_simd_kernels &&
253 nbnxn_simd_supported(mdlog, ir))
255 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl, hardwareInfo);
259 *kernel_type = nbnxnk4x4_PlainC;
265 GMX_LOG(mdlog.info).asParagraph().appendTextFormatted(
266 "Using %s %dx%d nonbonded short-range kernels",
267 lookup_kernel_name(*kernel_type),
268 nbnxn_kernel_to_cluster_i_size(*kernel_type),
269 nbnxn_kernel_to_cluster_j_size(*kernel_type));
271 if (nbnxnk4x4_PlainC == *kernel_type ||
272 nbnxnk8x8x8_PlainC == *kernel_type)
274 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
275 "WARNING: Using the slow %s kernels. This should\n"
276 "not happen during routine usage on supported platforms.",
277 lookup_kernel_name(*kernel_type));
282 void init_nb_verlet(const gmx::MDLogger &mdlog,
283 nonbonded_verlet_t **nb_verlet,
284 gmx_bool bFEP_NonBonded,
285 const t_inputrec *ir,
286 const t_forcerec *fr,
288 const gmx_hw_info_t &hardwareInfo,
289 const gmx_device_info_t *deviceInfo,
290 const gmx_mtop_t *mtop,
293 nonbonded_verlet_t *nbv = new nonbonded_verlet_t();
295 const EmulateGpuNonbonded emulateGpu =
296 ((getenv("GMX_EMULATE_GPU") != nullptr) ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
297 bool useGpu = deviceInfo != nullptr;
299 GMX_RELEASE_ASSERT(!(emulateGpu == EmulateGpuNonbonded::Yes && useGpu), "When GPU emulation is active, there cannot be a GPU assignment");
303 pick_nbnxn_kernel(mdlog, fr->use_simd_kernels, hardwareInfo,
304 useGpu, emulateGpu, ir,
306 &nbv->ewaldExclusionType_,
309 const bool haveMultipleDomains = (DOMAINDECOMP(cr) && cr->dd->nnodes > 1);
311 const bool pairlistIsSimple = nbv->pairlistIsSimple();
312 for (nbnxn_pairlist_set_t &pairlistSet : nbv->pairlistSets)
314 // TODO Change this to a constructor
315 /* The second parameter tells whether lists should be combined,
316 * this is currently only and always done for GPU lists.
318 nbnxn_init_pairlist_set(&pairlistSet, pairlistIsSimple, !pairlistIsSimple);
321 nbv->min_ci_balanced = 0;
323 nbv->listParams = std::make_unique<NbnxnListParameters>(ir->rlist);
324 setupDynamicPairlistPruning(mdlog, ir, mtop, box, nbv->kernelType_, fr->ic,
325 nbv->listParams.get());
327 nbv->nbs = std::make_unique<nbnxn_search>(ir->ePBC,
328 DOMAINDECOMP(cr) ? &cr->dd->nc : nullptr,
329 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : nullptr,
331 gmx_omp_nthreads_get(emntPairsearch));
333 int enbnxninitcombrule;
334 if (fr->ic->vdwtype == evdwCUT &&
335 (fr->ic->vdw_modifier == eintmodNONE ||
336 fr->ic->vdw_modifier == eintmodPOTSHIFT) &&
337 getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
339 /* Plain LJ cut-off: we can optimize with combination rules */
340 enbnxninitcombrule = enbnxninitcombruleDETECT;
342 else if (fr->ic->vdwtype == evdwPME)
344 /* LJ-PME: we need to use a combination rule for the grid */
345 if (fr->ljpme_combination_rule == eljpmeGEOM)
347 enbnxninitcombrule = enbnxninitcombruleGEOM;
351 enbnxninitcombrule = enbnxninitcombruleLB;
356 /* We use a full combination matrix: no rule required */
357 enbnxninitcombrule = enbnxninitcombruleNONE;
360 nbv->nbat = new nbnxn_atomdata_t(useGpu ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned);
361 int mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
362 if (ir->opts.ngener - ir->nwall == 1)
364 /* We have only one non-wall energy group, we do not need energy group
365 * support in the non-bondeds kernels, since all non-bonded energy
366 * contributions go to the first element of the energy group matrix.
368 mimimumNumEnergyGroupNonbonded = 1;
370 nbnxn_atomdata_init(mdlog,
375 mimimumNumEnergyGroupNonbonded,
376 pairlistIsSimple ? gmx_omp_nthreads_get(emntNonbonded) : 1);
380 /* init the NxN GPU data; the last argument tells whether we'll have
381 * both local and non-local NB calculation on GPU */
382 gpu_init(&nbv->gpu_nbv,
385 nbv->listParams.get(),
388 haveMultipleDomains);
390 if (const char *env = getenv("GMX_NB_MIN_CI"))
394 nbv->min_ci_balanced = strtol(env, &end, 10);
395 if (!end || (*end != 0) || nbv->min_ci_balanced < 0)
397 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
402 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
403 nbv->min_ci_balanced);
408 nbv->min_ci_balanced = gpu_min_ci_balanced(nbv->gpu_nbv);
411 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
412 nbv->min_ci_balanced);