2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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/pairlist_tuning.h"
56 #include "gromacs/simd/simd.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/logger.h"
61 #include "nbnxm_geometry.h"
62 #include "nbnxm_simd.h"
64 #include "pairlistset.h"
65 #include "pairlistsets.h"
66 #include "pairsearch.h"
71 /*! \brief Resources that can be used to execute non-bonded kernels on */
72 enum class NonbondedResource : int
79 /*! \brief Returns whether CPU SIMD support exists for the given inputrec
81 * If the return value is FALSE and fplog/cr != NULL, prints a fallback
82 * message to fplog/stderr.
84 static gmx_bool nbnxn_simd_supported(const gmx::MDLogger& mdlog, const t_inputrec* ir)
86 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
88 /* LJ PME with LB combination rule does 7 mesh operations.
89 * This so slow that we don't compile SIMD non-bonded kernels
91 GMX_LOG(mdlog.warning)
94 "LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling "
95 "back to plain C kernels");
102 /*! \brief Returns the most suitable CPU kernel type and Ewald handling */
103 static KernelSetup pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused* ir,
104 const gmx_hw_info_t gmx_unused& hardwareInfo)
106 KernelSetup kernelSetup;
110 kernelSetup.kernelType = KernelType::Cpu4x4_PlainC;
111 kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
115 #ifdef GMX_NBNXN_SIMD_4XN
116 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_4xN;
118 #ifdef GMX_NBNXN_SIMD_2XNN
119 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
122 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
123 /* We need to choose if we want 2x(N+N) or 4xN kernels.
124 * This is based on the SIMD acceleration choice and CPU information
125 * detected at runtime.
127 * 4xN calculates more (zero) interactions, but has less pair-search
128 * work and much better kernel instruction scheduling.
130 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
131 * which doesn't have FMA, both the analytical and tabulated Ewald
132 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
133 * 2x(4+4) because it results in significantly fewer pairs.
134 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
135 * 10% with HT, 50% without HT. As we currently don't detect the actual
136 * use of HT, use 4x8 to avoid a potential performance hit.
137 * On Intel Haswell 4x8 is always faster.
139 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_4xN;
141 if (!GMX_SIMD_HAVE_FMA && (EEL_PME_EWALD(ir->coulombtype) || EVDW_PME(ir->vdwtype)))
143 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
144 * There are enough instructions to make 2x(4+4) efficient.
146 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
149 if (hardwareInfo.haveAmdZen1Cpu)
151 /* One 256-bit FMA per cycle makes 2xNN faster */
152 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
154 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
157 if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
159 #ifdef GMX_NBNXN_SIMD_4XN
160 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_4xN;
163 "SIMD 4xN kernels requested, but GROMACS has been compiled without support "
164 "for these kernels");
167 if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
169 #ifdef GMX_NBNXN_SIMD_2XNN
170 kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
173 "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without "
174 "support for these kernels");
178 /* Analytical Ewald exclusion correction is only an option in
180 * Since table lookup's don't parallelize with SIMD, analytical
181 * will probably always be faster for a SIMD width of 8 or more.
182 * With FMA analytical is sometimes faster for a width if 4 as well.
183 * In single precision, this is faster on Bulldozer.
184 * On AMD Zen, tabulated Ewald kernels are faster on all 4 combinations
185 * of single or double precision and 128 or 256-bit AVX2.
189 (GMX_SIMD_REAL_WIDTH >= 8 || (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)) &&
191 !hardwareInfo.haveAmdZen1Cpu)
193 kernelSetup.ewaldExclusionType = EwaldExclusionType::Analytical;
197 kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
199 if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
201 kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
203 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
205 kernelSetup.ewaldExclusionType = EwaldExclusionType::Analytical;
212 const char* lookup_kernel_name(const KernelType kernelType)
214 const char* returnvalue = nullptr;
217 case KernelType::NotSet: returnvalue = "not set"; break;
218 case KernelType::Cpu4x4_PlainC: returnvalue = "plain C"; break;
219 case KernelType::Cpu4xN_Simd_4xN:
220 case KernelType::Cpu4xN_Simd_2xNN:
222 returnvalue = "SIMD";
224 returnvalue = "not available";
227 case KernelType::Gpu8x8x8: returnvalue = "GPU"; break;
228 case KernelType::Cpu8x8x8_PlainC: returnvalue = "plain C"; break;
230 default: gmx_fatal(FARGS, "Illegal kernel type selected");
235 /*! \brief Returns the most suitable kernel type and Ewald handling */
236 static KernelSetup pick_nbnxn_kernel(const gmx::MDLogger& mdlog,
237 gmx_bool use_simd_kernels,
238 const gmx_hw_info_t& hardwareInfo,
239 const NonbondedResource& nonbondedResource,
240 const t_inputrec* ir,
241 gmx_bool bDoNonbonded)
243 KernelSetup kernelSetup;
245 if (nonbondedResource == NonbondedResource::EmulateGpu)
247 kernelSetup.kernelType = KernelType::Cpu8x8x8_PlainC;
248 kernelSetup.ewaldExclusionType = EwaldExclusionType::DecidedByGpuModule;
252 GMX_LOG(mdlog.warning)
254 .appendText("Emulating a GPU run on the CPU (slow)");
257 else if (nonbondedResource == NonbondedResource::Gpu)
259 kernelSetup.kernelType = KernelType::Gpu8x8x8;
260 kernelSetup.ewaldExclusionType = EwaldExclusionType::DecidedByGpuModule;
264 if (use_simd_kernels && nbnxn_simd_supported(mdlog, ir))
266 kernelSetup = pick_nbnxn_kernel_cpu(ir, hardwareInfo);
270 kernelSetup.kernelType = KernelType::Cpu4x4_PlainC;
271 kernelSetup.ewaldExclusionType = EwaldExclusionType::Analytical;
279 .appendTextFormatted("Using %s %dx%d nonbonded short-range kernels",
280 lookup_kernel_name(kernelSetup.kernelType),
281 IClusterSizePerKernelType[kernelSetup.kernelType],
282 JClusterSizePerKernelType[kernelSetup.kernelType]);
284 if (KernelType::Cpu4x4_PlainC == kernelSetup.kernelType
285 || KernelType::Cpu8x8x8_PlainC == kernelSetup.kernelType)
287 GMX_LOG(mdlog.warning)
289 .appendTextFormatted(
290 "WARNING: Using the slow %s kernels. This should\n"
291 "not happen during routine usage on supported platforms.",
292 lookup_kernel_name(kernelSetup.kernelType));
296 GMX_RELEASE_ASSERT(kernelSetup.kernelType != KernelType::NotSet
297 && kernelSetup.ewaldExclusionType != EwaldExclusionType::NotSet,
298 "All kernel setup parameters should be set here");
305 PairlistSets::PairlistSets(const PairlistParams& pairlistParams,
306 const bool haveMultipleDomains,
307 const int minimumIlistCountForGpuBalancing) :
308 params_(pairlistParams),
309 minimumIlistCountForGpuBalancing_(minimumIlistCountForGpuBalancing)
311 localSet_ = std::make_unique<PairlistSet>(gmx::InteractionLocality::Local, params_);
313 if (haveMultipleDomains)
315 nonlocalSet_ = std::make_unique<PairlistSet>(gmx::InteractionLocality::NonLocal, params_);
322 /*! \brief Gets and returns the minimum i-list count for balacing based on the GPU used or env.var. when set */
323 static int getMinimumIlistCountForGpuBalancing(gmx_nbnxm_gpu_t* nbnxmGpu)
325 int minimumIlistCount;
327 if (const char* env = getenv("GMX_NB_MIN_CI"))
331 minimumIlistCount = strtol(env, &end, 10);
332 if (!end || (*end != 0) || minimumIlistCount < 0)
335 "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
340 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
346 minimumIlistCount = gpu_min_ci_balanced(nbnxmGpu);
350 "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU "
351 "multi-processors)\n",
356 return minimumIlistCount;
359 std::unique_ptr<nonbonded_verlet_t> init_nb_verlet(const gmx::MDLogger& mdlog,
360 gmx_bool bFEP_NonBonded,
361 const t_inputrec* ir,
362 const t_forcerec* fr,
364 const gmx_hw_info_t& hardwareInfo,
365 const gmx_device_info_t* deviceInfo,
366 const gmx_mtop_t* mtop,
368 gmx_wallcycle* wcycle)
370 const bool emulateGpu = (getenv("GMX_EMULATE_GPU") != nullptr);
371 const bool useGpu = deviceInfo != nullptr;
373 GMX_RELEASE_ASSERT(!(emulateGpu && useGpu),
374 "When GPU emulation is active, there cannot be a GPU assignment");
376 NonbondedResource nonbondedResource;
379 nonbondedResource = NonbondedResource::Gpu;
383 nonbondedResource = NonbondedResource::EmulateGpu;
387 nonbondedResource = NonbondedResource::Cpu;
390 Nbnxm::KernelSetup kernelSetup = pick_nbnxn_kernel(mdlog, fr->use_simd_kernels, hardwareInfo,
391 nonbondedResource, ir, fr->bNonbonded);
393 const bool haveMultipleDomains = (DOMAINDECOMP(cr) && cr->dd->nnodes > 1);
395 PairlistParams pairlistParams(kernelSetup.kernelType, bFEP_NonBonded, ir->rlist,
396 havePPDomainDecomposition(cr));
398 setupDynamicPairlistPruning(mdlog, ir, mtop, box, fr->ic, &pairlistParams);
400 int enbnxninitcombrule;
401 if (fr->ic->vdwtype == evdwCUT
402 && (fr->ic->vdw_modifier == eintmodNONE || fr->ic->vdw_modifier == eintmodPOTSHIFT)
403 && getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
405 /* Plain LJ cut-off: we can optimize with combination rules */
406 enbnxninitcombrule = enbnxninitcombruleDETECT;
408 else if (fr->ic->vdwtype == evdwPME)
410 /* LJ-PME: we need to use a combination rule for the grid */
411 if (fr->ljpme_combination_rule == eljpmeGEOM)
413 enbnxninitcombrule = enbnxninitcombruleGEOM;
417 enbnxninitcombrule = enbnxninitcombruleLB;
422 /* We use a full combination matrix: no rule required */
423 enbnxninitcombrule = enbnxninitcombruleNONE;
426 auto pinPolicy = (useGpu ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned);
428 auto nbat = std::make_unique<nbnxn_atomdata_t>(pinPolicy);
430 int mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
431 if (ir->opts.ngener - ir->nwall == 1)
433 /* We have only one non-wall energy group, we do not need energy group
434 * support in the non-bondeds kernels, since all non-bonded energy
435 * contributions go to the first element of the energy group matrix.
437 mimimumNumEnergyGroupNonbonded = 1;
439 nbnxn_atomdata_init(mdlog, nbat.get(), kernelSetup.kernelType, enbnxninitcombrule, fr->ntype,
440 fr->nbfp, mimimumNumEnergyGroupNonbonded,
441 (useGpu || emulateGpu) ? 1 : gmx_omp_nthreads_get(emntNonbonded));
443 gmx_nbnxm_gpu_t* gpu_nbv = nullptr;
444 int minimumIlistCountForGpuBalancing = 0;
447 /* init the NxN GPU data; the last argument tells whether we'll have
448 * both local and non-local NB calculation on GPU */
449 gpu_nbv = gpu_init(deviceInfo, fr->ic, pairlistParams, nbat.get(), cr->nodeid, haveMultipleDomains);
451 minimumIlistCountForGpuBalancing = getMinimumIlistCountForGpuBalancing(gpu_nbv);
454 auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, haveMultipleDomains,
455 minimumIlistCountForGpuBalancing);
457 auto pairSearch = std::make_unique<PairSearch>(
458 ir->pbcType, EI_TPI(ir->eI), DOMAINDECOMP(cr) ? &cr->dd->numCells : nullptr,
459 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : nullptr, pairlistParams.pairlistType,
460 bFEP_NonBonded, gmx_omp_nthreads_get(emntPairsearch), pinPolicy);
462 return std::make_unique<nonbonded_verlet_t>(std::move(pairlistSets), std::move(pairSearch),
463 std::move(nbat), kernelSetup, gpu_nbv, wcycle);
468 nonbonded_verlet_t::nonbonded_verlet_t(std::unique_ptr<PairlistSets> pairlistSets,
469 std::unique_ptr<PairSearch> pairSearch,
470 std::unique_ptr<nbnxn_atomdata_t> nbat_in,
471 const Nbnxm::KernelSetup& kernelSetup,
472 gmx_nbnxm_gpu_t* gpu_nbv_ptr,
473 gmx_wallcycle* wcycle) :
474 pairlistSets_(std::move(pairlistSets)),
475 pairSearch_(std::move(pairSearch)),
476 nbat(std::move(nbat_in)),
477 kernelSetup_(kernelSetup),
481 GMX_RELEASE_ASSERT(pairlistSets_, "Need valid pairlistSets");
482 GMX_RELEASE_ASSERT(pairSearch_, "Need valid search object");
483 GMX_RELEASE_ASSERT(nbat, "Need valid atomdata object");
486 nonbonded_verlet_t::~nonbonded_verlet_t()
488 Nbnxm::gpu_free(gpu_nbv);