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"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/mdlib/dispersioncorrection.h"
52 #include "gromacs/mdlib/force_flags.h"
53 #include "gromacs/mdlib/forcerec.h"
54 #include "gromacs/mdlib/gmx_omp_nthreads.h"
55 #include "gromacs/mdtypes/enerdata.h"
56 #include "gromacs/mdtypes/forcerec.h"
57 #include "gromacs/mdtypes/interaction_const.h"
58 #include "gromacs/mdtypes/mdatom.h"
59 #include "gromacs/mdtypes/simulation_workload.h"
60 #include "gromacs/nbnxm/atomdata.h"
61 #include "gromacs/nbnxm/gridset.h"
62 #include "gromacs/nbnxm/nbnxm.h"
63 #include "gromacs/nbnxm/nbnxm_simd.h"
64 #include "gromacs/nbnxm/pairlistset.h"
65 #include "gromacs/nbnxm/pairlistsets.h"
66 #include "gromacs/nbnxm/pairsearch.h"
67 #include "gromacs/pbcutil/ishift.h"
68 #include "gromacs/pbcutil/pbc.h"
69 #include "gromacs/simd/simd.h"
70 #include "gromacs/timing/cyclecounter.h"
71 #include "gromacs/utility/enumerationhelpers.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/logger.h"
75 #include "bench_system.h"
80 /*! \brief Checks the kernel setup
82 * Returns an error string when the kernel is not available.
84 static std::optional<std::string> checkKernelSetup(const KernelBenchOptions& options)
86 GMX_RELEASE_ASSERT(options.nbnxmSimd < BenchMarkKernels::Count
87 && options.nbnxmSimd != BenchMarkKernels::SimdAuto,
88 "Need a valid kernel SIMD type");
91 if ((options.nbnxmSimd != BenchMarkKernels::SimdNo && !GMX_SIMD)
92 #ifndef GMX_NBNXN_SIMD_4XN
93 || options.nbnxmSimd == BenchMarkKernels::Simd4XM
95 #ifndef GMX_NBNXN_SIMD_2XNN
96 || options.nbnxmSimd == BenchMarkKernels::Simd2XMM
100 return "the requested SIMD kernel was not set up at configuration time";
106 //! Helper to translate between the different enumeration values.
107 static KernelType translateBenchmarkEnum(const BenchMarkKernels& kernel)
109 int kernelInt = static_cast<int>(kernel);
110 return static_cast<KernelType>(kernelInt);
113 /*! \brief Returns the kernel setup
115 static KernelSetup getKernelSetup(const KernelBenchOptions& options)
117 auto messageWhenInvalid = checkKernelSetup(options);
118 GMX_RELEASE_ASSERT(!messageWhenInvalid, "Need valid options");
120 KernelSetup kernelSetup;
122 // The int enum options.nbnxnSimd is set up to match KernelType + 1
123 kernelSetup.kernelType = translateBenchmarkEnum(options.nbnxmSimd);
124 // The plain-C kernel does not support analytical ewald correction
125 if (kernelSetup.kernelType == KernelType::Cpu4x4_PlainC)
127 kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
131 kernelSetup.ewaldExclusionType = options.useTabulatedEwaldCorr ? EwaldExclusionType::Table
132 : EwaldExclusionType::Analytical;
138 //! Return an interaction constants struct with members used in the benchmark set appropriately
139 static interaction_const_t setupInteractionConst(const KernelBenchOptions& options)
142 interaction_const_t ic;
144 ic.vdwtype = evdwCUT;
145 ic.vdw_modifier = eintmodPOTSHIFT;
146 ic.rvdw = options.pairlistCutoff;
148 ic.eeltype = (options.coulombType == BenchMarkCoulomb::Pme ? eelPME : eelRF);
149 ic.coulomb_modifier = eintmodPOTSHIFT;
150 ic.rcoulomb = options.pairlistCutoff;
152 // Reaction-field with epsilon_rf=inf
153 // TODO: Replace by calc_rffac() after refactoring that
154 ic.k_rf = 0.5 * std::pow(ic.rcoulomb, -3);
155 ic.c_rf = 1 / ic.rcoulomb + ic.k_rf * ic.rcoulomb * ic.rcoulomb;
157 if (EEL_PME_EWALD(ic.eeltype))
159 // Ewald coefficients, we ignore the potential shift
160 GMX_RELEASE_ASSERT(options.ewaldcoeff_q > 0, "Ewald coefficient should be > 0");
161 ic.ewaldcoeff_q = options.ewaldcoeff_q;
162 ic.coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
163 init_interaction_const_tables(nullptr, &ic, 0);
169 //! Sets up and returns a Nbnxm object for the given benchmark options and system
170 static std::unique_ptr<nonbonded_verlet_t> setupNbnxmForBenchInstance(const KernelBenchOptions& options,
171 const gmx::BenchmarkSystem& system)
173 const auto pinPolicy = (options.useGpu ? gmx::PinningPolicy::PinnedIfSupported
174 : gmx::PinningPolicy::CannotBePinned);
175 const int numThreads = options.numThreads;
176 // Note: the options and Nbnxm combination rule enums values should match
177 const int combinationRule = static_cast<int>(options.ljCombinationRule);
179 auto messageWhenInvalid = checkKernelSetup(options);
180 if (messageWhenInvalid)
182 gmx_fatal(FARGS, "Requested kernel is unavailable because %s.", messageWhenInvalid->c_str());
184 Nbnxm::KernelSetup kernelSetup = getKernelSetup(options);
186 PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
188 GridSet gridSet(PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false,
189 numThreads, pinPolicy);
191 auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
194 std::make_unique<PairSearch>(PbcType::Xyz, false, nullptr, nullptr,
195 pairlistParams.pairlistType, false, numThreads, pinPolicy);
197 auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy);
199 // Put everything together
200 auto nbv = std::make_unique<nonbonded_verlet_t>(std::move(pairlistSets), std::move(pairSearch),
201 std::move(atomData), kernelSetup, nullptr, nullptr);
203 nbnxn_atomdata_init(gmx::MDLogger(), nbv->nbat.get(), kernelSetup.kernelType, combinationRule,
204 system.numAtomTypes, system.nonbondedParameters, 1, numThreads);
208 GMX_RELEASE_ASSERT(!TRICLINIC(system.box), "Only rectangular unit-cells are supported here");
209 const rvec lowerCorner = { 0, 0, 0 };
210 const rvec upperCorner = { system.box[XX][XX], system.box[YY][YY], system.box[ZZ][ZZ] };
212 gmx::ArrayRef<const int> atomInfo;
213 if (options.useHalfLJOptimization)
215 atomInfo = system.atomInfoOxygenVdw;
219 atomInfo = system.atomInfoAllVdw;
222 const real atomDensity = system.coordinates.size() / det(system.box);
224 nbnxn_put_on_grid(nbv.get(), system.box, 0, lowerCorner, upperCorner, nullptr,
225 { 0, int(system.coordinates.size()) }, atomDensity, atomInfo,
226 system.coordinates, 0, nullptr);
228 nbv->constructPairlist(gmx::InteractionLocality::Local, system.excls, 0, &nrnb);
230 nbv->setAtomProperties(system.atomTypes, system.charges, atomInfo);
235 //! Add the options instance to the list for all requested kernel SIMD types
236 static void expandSimdOptionAndPushBack(const KernelBenchOptions& options,
237 std::vector<KernelBenchOptions>* optionsList)
239 if (options.nbnxmSimd == BenchMarkKernels::SimdAuto)
241 bool addedInstance = false;
242 #ifdef GMX_NBNXN_SIMD_4XN
243 optionsList->push_back(options);
244 optionsList->back().nbnxmSimd = BenchMarkKernels::Simd4XM;
245 addedInstance = true;
247 #ifdef GMX_NBNXN_SIMD_2XNN
248 optionsList->push_back(options);
249 optionsList->back().nbnxmSimd = BenchMarkKernels::Simd2XMM;
250 addedInstance = true;
254 optionsList->push_back(options);
255 optionsList->back().nbnxmSimd = BenchMarkKernels::SimdNo;
260 optionsList->push_back(options);
264 //! Sets up and runs the requested benchmark instance and prints the results
266 // When \p doWarmup is true runs the warmup iterations instead
267 // of the normal ones and does not print any results
268 static void setupAndRunInstance(const gmx::BenchmarkSystem& system,
269 const KernelBenchOptions& options,
272 // Generate an, accurate, estimate of the number of non-zero pair interactions
273 const real atomDensity = system.coordinates.size() / det(system.box);
274 const real numPairsWithinCutoff =
275 atomDensity * 4.0 / 3.0 * M_PI * std::pow(options.pairlistCutoff, 3);
276 const real numUsefulPairs = system.coordinates.size() * 0.5 * (numPairsWithinCutoff + 1);
278 std::unique_ptr<nonbonded_verlet_t> nbv = setupNbnxmForBenchInstance(options, system);
280 // We set the interaction cut-off to the pairlist cut-off
281 interaction_const_t ic = setupInteractionConst(options);
285 gmx_enerdata_t enerd(1, 0);
287 gmx::StepWorkload stepWork;
288 stepWork.computeForces = true;
289 if (options.computeVirialAndEnergy)
291 stepWork.computeVirial = true;
292 stepWork.computeEnergy = true;
295 const gmx::EnumerationArray<BenchMarkKernels, std::string> kernelNames = { "auto", "no", "4xM",
298 const gmx::EnumerationArray<BenchMarkCombRule, std::string> combruleNames = { "geom.", "LB",
303 fprintf(stdout, "%-7s %-4s %-5s %-4s ",
304 options.coulombType == BenchMarkCoulomb::Pme ? "Ewald" : "RF",
305 options.useHalfLJOptimization ? "half" : "all",
306 combruleNames[options.ljCombinationRule].c_str(), kernelNames[options.nbnxmSimd].c_str());
309 // Run pre-iteration to avoid cache misses
310 for (int iter = 0; iter < options.numPreIterations; iter++)
312 nbv->dispatchNonbondedKernel(gmx::InteractionLocality::Local, ic, stepWork, enbvClearFYes,
313 system.forceRec, &enerd, &nrnb);
316 const int numIterations = (doWarmup ? options.numWarmupIterations : options.numIterations);
317 const PairlistSet& pairlistSet = nbv->pairlistSets().pairlistSet(gmx::InteractionLocality::Local);
318 const gmx::index numPairs = pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_ + pairlistSet.natpair_q_;
319 gmx_cycles_t cycles = gmx_cycles_read();
320 for (int iter = 0; iter < numIterations; iter++)
322 // Run the kernel without force clearing
323 nbv->dispatchNonbondedKernel(gmx::InteractionLocality::Local, ic, stepWork, enbvClearFNo,
324 system.forceRec, &enerd, &nrnb);
326 cycles = gmx_cycles_read() - cycles;
329 const double dCycles = static_cast<double>(cycles);
330 if (options.cyclesPerPair)
332 fprintf(stdout, "%10.3f %10.4f %8.4f %8.4f\n", cycles * 1e-6,
333 dCycles / options.numIterations * 1e-6, dCycles / (options.numIterations * numPairs),
334 dCycles / (options.numIterations * numUsefulPairs));
338 fprintf(stdout, "%10.3f %10.4f %8.4f %8.4f\n", dCycles * 1e-6,
339 dCycles / options.numIterations * 1e-6, options.numIterations * numPairs / dCycles,
340 options.numIterations * numUsefulPairs / dCycles);
345 void bench(const int sizeFactor, const KernelBenchOptions& options)
347 // We don't want to call gmx_omp_nthreads_init(), so we init what we need
348 gmx_omp_nthreads_set(emntPairsearch, options.numThreads);
349 gmx_omp_nthreads_set(emntNonbonded, options.numThreads);
351 const gmx::BenchmarkSystem system(sizeFactor);
353 real minBoxSize = norm(system.box[XX]);
354 for (int dim = YY; dim < DIM; dim++)
356 minBoxSize = std::min(minBoxSize, norm(system.box[dim]));
358 if (options.pairlistCutoff > 0.5 * minBoxSize)
360 gmx_fatal(FARGS, "The cut-off should be shorter than half the box size");
363 std::vector<KernelBenchOptions> optionsList;
366 KernelBenchOptions opt = options;
367 gmx::EnumerationWrapper<BenchMarkCoulomb> coulombIter;
368 for (auto coulombType : coulombIter)
370 opt.coulombType = coulombType;
371 for (int halfLJ = 0; halfLJ <= 1; halfLJ++)
373 opt.useHalfLJOptimization = (halfLJ == 1);
375 gmx::EnumerationWrapper<BenchMarkCombRule> combRuleIter;
376 for (auto combRule : combRuleIter)
378 if (combRule == BenchMarkCombRule::RuleNone)
382 opt.ljCombinationRule = combRule;
384 expandSimdOptionAndPushBack(opt, &optionsList);
391 expandSimdOptionAndPushBack(options, &optionsList);
393 GMX_RELEASE_ASSERT(!optionsList.empty(), "Expect at least on benchmark setup");
396 if (options.nbnxmSimd != BenchMarkKernels::SimdNo)
398 fprintf(stdout, "SIMD width: %d\n", GMX_SIMD_REAL_WIDTH);
401 fprintf(stdout, "System size: %zu atoms\n", system.coordinates.size());
402 fprintf(stdout, "Cut-off radius: %g nm\n", options.pairlistCutoff);
403 fprintf(stdout, "Number of threads: %d\n", options.numThreads);
404 fprintf(stdout, "Number of iterations: %d\n", options.numIterations);
405 fprintf(stdout, "Compute energies: %s\n", options.computeVirialAndEnergy ? "yes" : "no");
406 if (options.coulombType != BenchMarkCoulomb::ReactionField)
408 fprintf(stdout, "Ewald excl. corr.: %s\n",
409 options.nbnxmSimd == BenchMarkKernels::SimdNo || options.useTabulatedEwaldCorr
415 if (options.numWarmupIterations > 0)
417 setupAndRunInstance(system, optionsList[0], true);
420 fprintf(stdout, "Coulomb LJ comb. SIMD Mcycles Mcycles/it. %s\n",
421 options.cyclesPerPair ? "cycles/pair" : "pairs/cycle");
422 fprintf(stdout, " total useful\n");
424 for (const auto& optionsInstance : optionsList)
426 setupAndRunInstance(system, optionsInstance, false);