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.
38 * This file defines functions for setting up kernel benchmarks
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_nbnxm
46 #include "bench_setup.h"
48 #include "gromacs/compat/optional.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/mdlib/dispersioncorrection.h"
51 #include "gromacs/mdlib/force_flags.h"
52 #include "gromacs/mdlib/forcerec.h"
53 #include "gromacs/mdlib/gmx_omp_nthreads.h"
54 #include "gromacs/mdtypes/enerdata.h"
55 #include "gromacs/mdtypes/forcerec.h"
56 #include "gromacs/mdtypes/interaction_const.h"
57 #include "gromacs/mdtypes/mdatom.h"
58 #include "gromacs/mdtypes/simulation_workload.h"
59 #include "gromacs/nbnxm/atomdata.h"
60 #include "gromacs/nbnxm/gridset.h"
61 #include "gromacs/nbnxm/nbnxm.h"
62 #include "gromacs/nbnxm/nbnxm_simd.h"
63 #include "gromacs/nbnxm/pairlistset.h"
64 #include "gromacs/nbnxm/pairlistsets.h"
65 #include "gromacs/nbnxm/pairsearch.h"
66 #include "gromacs/pbcutil/ishift.h"
67 #include "gromacs/pbcutil/pbc.h"
68 #include "gromacs/simd/simd.h"
69 #include "gromacs/timing/cyclecounter.h"
70 #include "gromacs/utility/enumerationhelpers.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/logger.h"
74 #include "bench_system.h"
79 /*! \brief Checks the kernel setup
81 * Returns an error string when the kernel is not available.
83 static gmx::compat::optional<std::string> checkKernelSetup(const KernelBenchOptions& options)
85 GMX_RELEASE_ASSERT(options.nbnxmSimd < BenchMarkKernels::Count
86 && options.nbnxmSimd != BenchMarkKernels::SimdAuto,
87 "Need a valid kernel SIMD type");
90 if ((options.nbnxmSimd != BenchMarkKernels::SimdNo && !GMX_SIMD)
91 #ifndef GMX_NBNXN_SIMD_4XN
92 || options.nbnxmSimd == BenchMarkKernels::Simd4XM
94 #ifndef GMX_NBNXN_SIMD_2XNN
95 || options.nbnxmSimd == BenchMarkKernels::Simd2XMM
99 return "the requested SIMD kernel was not set up at configuration time";
105 //! Helper to translate between the different enumeration values.
106 static KernelType translateBenchmarkEnum(const BenchMarkKernels& kernel)
108 int kernelInt = static_cast<int>(kernel);
109 return static_cast<KernelType>(kernelInt);
112 /*! \brief Returns the kernel setup
114 static KernelSetup getKernelSetup(const KernelBenchOptions& options)
116 auto messageWhenInvalid = checkKernelSetup(options);
117 GMX_RELEASE_ASSERT(!messageWhenInvalid, "Need valid options");
119 KernelSetup kernelSetup;
121 // The int enum options.nbnxnSimd is set up to match KernelType + 1
122 kernelSetup.kernelType = translateBenchmarkEnum(options.nbnxmSimd);
123 // The plain-C kernel does not support analytical ewald correction
124 if (kernelSetup.kernelType == KernelType::Cpu4x4_PlainC)
126 kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
130 kernelSetup.ewaldExclusionType = options.useTabulatedEwaldCorr ? EwaldExclusionType::Table
131 : EwaldExclusionType::Analytical;
137 //! Return an interaction constants struct with members used in the benchmark set appropriately
138 static interaction_const_t setupInteractionConst(const KernelBenchOptions& options)
141 interaction_const_t ic;
143 ic.vdwtype = evdwCUT;
144 ic.vdw_modifier = eintmodPOTSHIFT;
145 ic.rvdw = options.pairlistCutoff;
147 ic.eeltype = (options.coulombType == BenchMarkCoulomb::Pme ? eelPME : eelRF);
148 ic.coulomb_modifier = eintmodPOTSHIFT;
149 ic.rcoulomb = options.pairlistCutoff;
151 // Reaction-field with epsilon_rf=inf
152 // TODO: Replace by calc_rffac() after refactoring that
153 ic.k_rf = 0.5 * std::pow(ic.rcoulomb, -3);
154 ic.c_rf = 1 / ic.rcoulomb + ic.k_rf * ic.rcoulomb * ic.rcoulomb;
156 if (EEL_PME_EWALD(ic.eeltype))
158 // Ewald coefficients, we ignore the potential shift
159 GMX_RELEASE_ASSERT(options.ewaldcoeff_q > 0, "Ewald coefficient should be > 0");
160 ic.ewaldcoeff_q = options.ewaldcoeff_q;
161 ic.coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
162 init_interaction_const_tables(nullptr, &ic);
168 //! Sets up and returns a Nbnxm object for the given benchmark options and system
169 static std::unique_ptr<nonbonded_verlet_t> setupNbnxmForBenchInstance(const KernelBenchOptions& options,
170 const gmx::BenchmarkSystem& system)
172 const auto pinPolicy = (options.useGpu ? gmx::PinningPolicy::PinnedIfSupported
173 : gmx::PinningPolicy::CannotBePinned);
174 const int numThreads = options.numThreads;
175 // Note: the options and Nbnxm combination rule enums values should match
176 const int combinationRule = static_cast<int>(options.ljCombinationRule);
178 auto messageWhenInvalid = checkKernelSetup(options);
179 if (messageWhenInvalid)
181 gmx_fatal(FARGS, "Requested kernel is unavailable because %s.", messageWhenInvalid->c_str());
183 Nbnxm::KernelSetup kernelSetup = getKernelSetup(options);
185 PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
187 GridSet gridSet(PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false,
188 numThreads, pinPolicy);
190 auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
193 std::make_unique<PairSearch>(PbcType::Xyz, false, nullptr, nullptr,
194 pairlistParams.pairlistType, false, numThreads, pinPolicy);
196 auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy);
198 // Put everything together
199 auto nbv = std::make_unique<nonbonded_verlet_t>(std::move(pairlistSets), std::move(pairSearch),
200 std::move(atomData), kernelSetup, nullptr, nullptr);
202 nbnxn_atomdata_init(gmx::MDLogger(), nbv->nbat.get(), kernelSetup.kernelType, combinationRule,
203 system.numAtomTypes, system.nonbondedParameters, 1, numThreads);
207 GMX_RELEASE_ASSERT(!TRICLINIC(system.box), "Only rectangular unit-cells are supported here");
208 const rvec lowerCorner = { 0, 0, 0 };
209 const rvec upperCorner = { system.box[XX][XX], system.box[YY][YY], system.box[ZZ][ZZ] };
211 gmx::ArrayRef<const int> atomInfo;
212 if (options.useHalfLJOptimization)
214 atomInfo = system.atomInfoOxygenVdw;
218 atomInfo = system.atomInfoAllVdw;
221 const real atomDensity = system.coordinates.size() / det(system.box);
223 nbnxn_put_on_grid(nbv.get(), system.box, 0, lowerCorner, upperCorner, nullptr,
224 { 0, int(system.coordinates.size()) }, atomDensity, atomInfo,
225 system.coordinates, 0, nullptr);
227 nbv->constructPairlist(gmx::InteractionLocality::Local, system.excls, 0, &nrnb);
230 // We only use (read) the atom type and charge from mdatoms
231 mdatoms.typeA = const_cast<int*>(system.atomTypes.data());
232 mdatoms.chargeA = const_cast<real*>(system.charges.data());
233 nbv->setAtomProperties(mdatoms, atomInfo);
238 //! Add the options instance to the list for all requested kernel SIMD types
239 static void expandSimdOptionAndPushBack(const KernelBenchOptions& options,
240 std::vector<KernelBenchOptions>* optionsList)
242 if (options.nbnxmSimd == BenchMarkKernels::SimdAuto)
244 bool addedInstance = false;
245 #ifdef GMX_NBNXN_SIMD_4XN
246 optionsList->push_back(options);
247 optionsList->back().nbnxmSimd = BenchMarkKernels::Simd4XM;
248 addedInstance = true;
250 #ifdef GMX_NBNXN_SIMD_2XNN
251 optionsList->push_back(options);
252 optionsList->back().nbnxmSimd = BenchMarkKernels::Simd2XMM;
253 addedInstance = true;
257 optionsList->push_back(options);
258 optionsList->back().nbnxmSimd = BenchMarkKernels::SimdNo;
263 optionsList->push_back(options);
267 //! Sets up and runs the requested benchmark instance and prints the results
269 // When \p doWarmup is true runs the warmup iterations instead
270 // of the normal ones and does not print any results
271 static void setupAndRunInstance(const gmx::BenchmarkSystem& system,
272 const KernelBenchOptions& options,
275 // Generate an, accurate, estimate of the number of non-zero pair interactions
276 const real atomDensity = system.coordinates.size() / det(system.box);
277 const real numPairsWithinCutoff =
278 atomDensity * 4.0 / 3.0 * M_PI * std::pow(options.pairlistCutoff, 3);
279 const real numUsefulPairs = system.coordinates.size() * 0.5 * (numPairsWithinCutoff + 1);
281 std::unique_ptr<nonbonded_verlet_t> nbv = setupNbnxmForBenchInstance(options, system);
283 // We set the interaction cut-off to the pairlist cut-off
284 interaction_const_t ic = setupInteractionConst(options);
288 gmx_enerdata_t enerd(1, 0);
290 gmx::StepWorkload stepWork;
291 stepWork.computeForces = true;
292 if (options.computeVirialAndEnergy)
294 stepWork.computeVirial = true;
295 stepWork.computeEnergy = true;
298 const gmx::EnumerationArray<BenchMarkKernels, std::string> kernelNames = { "auto", "no", "4xM",
301 const gmx::EnumerationArray<BenchMarkCombRule, std::string> combruleNames = { "geom.", "LB",
306 fprintf(stdout, "%-7s %-4s %-5s %-4s ",
307 options.coulombType == BenchMarkCoulomb::Pme ? "Ewald" : "RF",
308 options.useHalfLJOptimization ? "half" : "all",
309 combruleNames[options.ljCombinationRule].c_str(), kernelNames[options.nbnxmSimd].c_str());
312 // Run pre-iteration to avoid cache misses
313 for (int iter = 0; iter < options.numPreIterations; iter++)
315 nbv->dispatchNonbondedKernel(gmx::InteractionLocality::Local, ic, stepWork, enbvClearFYes,
316 system.forceRec, &enerd, &nrnb);
319 const int numIterations = (doWarmup ? options.numWarmupIterations : options.numIterations);
320 const PairlistSet& pairlistSet = nbv->pairlistSets().pairlistSet(gmx::InteractionLocality::Local);
321 const gmx::index numPairs = pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_ + pairlistSet.natpair_q_;
322 gmx_cycles_t cycles = gmx_cycles_read();
323 for (int iter = 0; iter < numIterations; iter++)
325 // Run the kernel without force clearing
326 nbv->dispatchNonbondedKernel(gmx::InteractionLocality::Local, ic, stepWork, enbvClearFNo,
327 system.forceRec, &enerd, &nrnb);
329 cycles = gmx_cycles_read() - cycles;
332 const double dCycles = static_cast<double>(cycles);
333 if (options.cyclesPerPair)
335 fprintf(stdout, "%10.3f %10.4f %8.4f %8.4f\n", cycles * 1e-6,
336 dCycles / options.numIterations * 1e-6, dCycles / (options.numIterations * numPairs),
337 dCycles / (options.numIterations * numUsefulPairs));
341 fprintf(stdout, "%10.3f %10.4f %8.4f %8.4f\n", dCycles * 1e-6,
342 dCycles / options.numIterations * 1e-6, options.numIterations * numPairs / dCycles,
343 options.numIterations * numUsefulPairs / dCycles);
348 void bench(const int sizeFactor, const KernelBenchOptions& options)
350 // We don't want to call gmx_omp_nthreads_init(), so we init what we need
351 gmx_omp_nthreads_set(emntPairsearch, options.numThreads);
352 gmx_omp_nthreads_set(emntNonbonded, options.numThreads);
354 const gmx::BenchmarkSystem system(sizeFactor);
356 real minBoxSize = norm(system.box[XX]);
357 for (int dim = YY; dim < DIM; dim++)
359 minBoxSize = std::min(minBoxSize, norm(system.box[dim]));
361 if (options.pairlistCutoff > 0.5 * minBoxSize)
363 gmx_fatal(FARGS, "The cut-off should be shorter than half the box size");
366 std::vector<KernelBenchOptions> optionsList;
369 KernelBenchOptions opt = options;
370 gmx::EnumerationWrapper<BenchMarkCoulomb> coulombIter;
371 for (auto coulombType : coulombIter)
373 opt.coulombType = coulombType;
374 for (int halfLJ = 0; halfLJ <= 1; halfLJ++)
376 opt.useHalfLJOptimization = (halfLJ == 1);
378 gmx::EnumerationWrapper<BenchMarkCombRule> combRuleIter;
379 for (auto combRule : combRuleIter)
381 if (combRule == BenchMarkCombRule::RuleNone)
385 opt.ljCombinationRule = combRule;
387 expandSimdOptionAndPushBack(opt, &optionsList);
394 expandSimdOptionAndPushBack(options, &optionsList);
396 GMX_RELEASE_ASSERT(!optionsList.empty(), "Expect at least on benchmark setup");
399 if (options.nbnxmSimd != BenchMarkKernels::SimdNo)
401 fprintf(stdout, "SIMD width: %d\n", GMX_SIMD_REAL_WIDTH);
404 fprintf(stdout, "System size: %zu atoms\n", system.coordinates.size());
405 fprintf(stdout, "Cut-off radius: %g nm\n", options.pairlistCutoff);
406 fprintf(stdout, "Number of threads: %d\n", options.numThreads);
407 fprintf(stdout, "Number of iterations: %d\n", options.numIterations);
408 fprintf(stdout, "Compute energies: %s\n", options.computeVirialAndEnergy ? "yes" : "no");
409 if (options.coulombType != BenchMarkCoulomb::ReactionField)
411 fprintf(stdout, "Ewald excl. corr.: %s\n",
412 options.nbnxmSimd == BenchMarkKernels::SimdNo || options.useTabulatedEwaldCorr
418 if (options.numWarmupIterations > 0)
420 setupAndRunInstance(system, optionsList[0], true);
423 fprintf(stdout, "Coulomb LJ comb. SIMD Mcycles Mcycles/it. %s\n",
424 options.cyclesPerPair ? "cycles/pair" : "pairs/cycle");
425 fprintf(stdout, " total useful\n");
427 for (const auto& optionsInstance : optionsList)
429 setupAndRunInstance(system, optionsInstance, false);