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.h"
58 #include "gromacs/nbnxm/pairlist_tuning.h"
59 #include "gromacs/nbnxm/pairlistset.h"
60 #include "gromacs/simd/simd.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/logger.h"
70 /*! \brief Returns whether CPU SIMD support exists for the given inputrec
72 * If the return value is FALSE and fplog/cr != NULL, prints a fallback
73 * message to fplog/stderr.
75 static gmx_bool nbnxn_simd_supported(const gmx::MDLogger &mdlog,
78 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
80 /* LJ PME with LB combination rule does 7 mesh operations.
81 * This so slow that we don't compile SIMD non-bonded kernels
83 GMX_LOG(mdlog.warning).asParagraph().appendText("LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels");
90 /*! \brief Returns the most suitable CPU kernel type and Ewald handling */
92 pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
93 const gmx_hw_info_t gmx_unused &hardwareInfo)
95 KernelSetup kernelSetup;
99 kernelSetup.kernelType = KernelType::Cpu4x4_PlainC;
100 kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
104 #ifdef GMX_NBNXN_SIMD_4XN
105 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_4xN;
107 #ifdef GMX_NBNXN_SIMD_2XNN
108 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
111 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
112 /* We need to choose if we want 2x(N+N) or 4xN kernels.
113 * This is based on the SIMD acceleration choice and CPU information
114 * detected at runtime.
116 * 4xN calculates more (zero) interactions, but has less pair-search
117 * work and much better kernel instruction scheduling.
119 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
120 * which doesn't have FMA, both the analytical and tabulated Ewald
121 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
122 * 2x(4+4) because it results in significantly fewer pairs.
123 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
124 * 10% with HT, 50% without HT. As we currently don't detect the actual
125 * use of HT, use 4x8 to avoid a potential performance hit.
126 * On Intel Haswell 4x8 is always faster.
128 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_4xN;
130 if (!GMX_SIMD_HAVE_FMA && (EEL_PME_EWALD(ir->coulombtype) ||
131 EVDW_PME(ir->vdwtype)))
133 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
134 * There are enough instructions to make 2x(4+4) efficient.
136 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
139 if (hardwareInfo.haveAmdZenCpu)
141 /* One 256-bit FMA per cycle makes 2xNN faster */
142 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
144 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
147 if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
149 #ifdef GMX_NBNXN_SIMD_4XN
150 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_4xN;
152 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
155 if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
157 #ifdef GMX_NBNXN_SIMD_2XNN
158 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
160 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
164 /* Analytical Ewald exclusion correction is only an option in
166 * Since table lookup's don't parallelize with SIMD, analytical
167 * will probably always be faster for a SIMD width of 8 or more.
168 * With FMA analytical is sometimes faster for a width if 4 as well.
169 * In single precision, this is faster on Bulldozer.
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.
175 (GMX_SIMD_REAL_WIDTH >= 8 ||
176 (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)) &&
178 !hardwareInfo.haveAmdZenCpu)
180 kernelSetup.ewaldExclusionType = EwaldExclusionType::Analytical;
184 kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
186 if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
188 kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
190 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
192 kernelSetup.ewaldExclusionType = EwaldExclusionType::Analytical;
200 const char *lookup_kernel_name(const KernelType kernelType)
202 const char *returnvalue = nullptr;
205 case KernelType::NotSet:
206 returnvalue = "not set";
208 case KernelType::Cpu4x4_PlainC:
209 returnvalue = "plain C";
211 case KernelType::Cpu4xN_Simd_4xN:
212 case KernelType::Cpu4xN_Simd_2xNN:
214 returnvalue = "SIMD";
216 returnvalue = "not available";
219 case KernelType::Gpu8x8x8: returnvalue = "GPU"; break;
220 case KernelType::Cpu8x8x8_PlainC: returnvalue = "plain C"; break;
223 gmx_fatal(FARGS, "Illegal kernel type selected");
228 /*! \brief Returns the most suitable kernel type and Ewald handling */
230 pick_nbnxn_kernel(const gmx::MDLogger &mdlog,
231 gmx_bool use_simd_kernels,
232 const gmx_hw_info_t &hardwareInfo,
233 const NonbondedResource &nonbondedResource,
234 const t_inputrec *ir,
235 gmx_bool bDoNonbonded)
237 KernelSetup kernelSetup;
239 if (nonbondedResource == NonbondedResource::EmulateGpu)
241 kernelSetup.kernelType = KernelType::Cpu8x8x8_PlainC;
242 kernelSetup.ewaldExclusionType = EwaldExclusionType::DecidedByGpuModule;
246 GMX_LOG(mdlog.warning).asParagraph().appendText("Emulating a GPU run on the CPU (slow)");
249 else if (nonbondedResource == NonbondedResource::Gpu)
251 kernelSetup.kernelType = KernelType::Gpu8x8x8;
252 kernelSetup.ewaldExclusionType = EwaldExclusionType::DecidedByGpuModule;
256 if (use_simd_kernels &&
257 nbnxn_simd_supported(mdlog, ir))
259 kernelSetup = pick_nbnxn_kernel_cpu(ir, hardwareInfo);
263 kernelSetup.kernelType = KernelType::Cpu4x4_PlainC;
264 kernelSetup.ewaldExclusionType = EwaldExclusionType::Analytical;
270 GMX_LOG(mdlog.info).asParagraph().appendTextFormatted(
271 "Using %s %dx%d nonbonded short-range kernels",
272 lookup_kernel_name(kernelSetup.kernelType),
273 IClusterSizePerKernelType[kernelSetup.kernelType],
274 JClusterSizePerKernelType[kernelSetup.kernelType]);
276 if (KernelType::Cpu4x4_PlainC == kernelSetup.kernelType ||
277 KernelType::Cpu8x8x8_PlainC == kernelSetup.kernelType)
279 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
280 "WARNING: Using the slow %s kernels. This should\n"
281 "not happen during routine usage on supported platforms.",
282 lookup_kernel_name(kernelSetup.kernelType));
286 GMX_RELEASE_ASSERT(kernelSetup.kernelType != KernelType::NotSet &&
287 kernelSetup.ewaldExclusionType != EwaldExclusionType::NotSet,
288 "All kernel setup parameters should be set here");
295 nonbonded_verlet_t::PairlistSets::PairlistSets(const NbnxnListParameters &listParams,
296 const bool haveMultipleDomains,
297 const int minimumIlistCountForGpuBalancing) :
299 minimumIlistCountForGpuBalancing_(minimumIlistCountForGpuBalancing)
301 localSet_ = std::make_unique<nbnxn_pairlist_set_t>(params_);
303 if (haveMultipleDomains)
305 nonlocalSet_ = std::make_unique<nbnxn_pairlist_set_t>(params_);
312 /*! \brief Gets and returns the minimum i-list count for balacing based on the GPU used or env.var. when set */
313 static int getMinimumIlistCountForGpuBalancing(gmx_nbnxn_gpu_t *nbnxmGpu)
315 int minimumIlistCount;
317 if (const char *env = getenv("GMX_NB_MIN_CI"))
321 minimumIlistCount = strtol(env, &end, 10);
322 if (!end || (*end != 0) || minimumIlistCount < 0)
324 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
329 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
335 minimumIlistCount = gpu_min_ci_balanced(nbnxmGpu);
338 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
343 return minimumIlistCount;
346 void init_nb_verlet(const gmx::MDLogger &mdlog,
347 nonbonded_verlet_t **nb_verlet,
348 gmx_bool bFEP_NonBonded,
349 const t_inputrec *ir,
350 const t_forcerec *fr,
352 const gmx_hw_info_t &hardwareInfo,
353 const gmx_device_info_t *deviceInfo,
354 const gmx_mtop_t *mtop,
357 nonbonded_verlet_t *nbv = new nonbonded_verlet_t();
359 const bool emulateGpu = (getenv("GMX_EMULATE_GPU") != nullptr);
360 const bool useGpu = deviceInfo != nullptr;
362 GMX_RELEASE_ASSERT(!(emulateGpu && useGpu), "When GPU emulation is active, there cannot be a GPU assignment");
364 NonbondedResource nonbondedResource;
367 nonbondedResource = NonbondedResource::Gpu;
371 nonbondedResource = NonbondedResource::EmulateGpu;
375 nonbondedResource = NonbondedResource::Cpu;
380 nbv->setKernelSetup(pick_nbnxn_kernel(mdlog, fr->use_simd_kernels, hardwareInfo,
381 nonbondedResource, ir,
384 const bool haveMultipleDomains = (DOMAINDECOMP(cr) && cr->dd->nnodes > 1);
386 NbnxnListParameters listParams(nbv->kernelSetup().kernelType, ir->rlist);
388 setupDynamicPairlistPruning(mdlog, ir, mtop, box, fr->ic,
391 int enbnxninitcombrule;
392 if (fr->ic->vdwtype == evdwCUT &&
393 (fr->ic->vdw_modifier == eintmodNONE ||
394 fr->ic->vdw_modifier == eintmodPOTSHIFT) &&
395 getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
397 /* Plain LJ cut-off: we can optimize with combination rules */
398 enbnxninitcombrule = enbnxninitcombruleDETECT;
400 else if (fr->ic->vdwtype == evdwPME)
402 /* LJ-PME: we need to use a combination rule for the grid */
403 if (fr->ljpme_combination_rule == eljpmeGEOM)
405 enbnxninitcombrule = enbnxninitcombruleGEOM;
409 enbnxninitcombrule = enbnxninitcombruleLB;
414 /* We use a full combination matrix: no rule required */
415 enbnxninitcombrule = enbnxninitcombruleNONE;
418 nbv->nbat = new nbnxn_atomdata_t(useGpu ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned);
419 int mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
420 if (ir->opts.ngener - ir->nwall == 1)
422 /* We have only one non-wall energy group, we do not need energy group
423 * support in the non-bondeds kernels, since all non-bonded energy
424 * contributions go to the first element of the energy group matrix.
426 mimimumNumEnergyGroupNonbonded = 1;
428 nbnxn_atomdata_init(mdlog,
430 nbv->kernelSetup().kernelType,
433 mimimumNumEnergyGroupNonbonded,
434 nbv->pairlistIsSimple() ? gmx_omp_nthreads_get(emntNonbonded) : 1);
436 int minimumIlistCountForGpuBalancing = 0;
439 /* init the NxN GPU data; the last argument tells whether we'll have
440 * both local and non-local NB calculation on GPU */
441 gpu_init(&nbv->gpu_nbv,
447 haveMultipleDomains);
449 minimumIlistCountForGpuBalancing = getMinimumIlistCountForGpuBalancing(nbv->gpu_nbv);
453 std::make_unique<nonbonded_verlet_t::PairlistSets>(listParams,
455 minimumIlistCountForGpuBalancing);
457 nbv->nbs = std::make_unique<nbnxn_search>(ir->ePBC,
458 DOMAINDECOMP(cr) ? &cr->dd->nc : nullptr,
459 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : nullptr,
461 gmx_omp_nthreads_get(emntPairsearch));