Refactor md_enums
[alexxy/gromacs.git] / src / gromacs / nbnxm / benchmark / bench_setup.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020,2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  * \brief
38  * This file defines functions for setting up kernel benchmarks
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \ingroup module_nbnxm
42  */
43
44 #include "gmxpre.h"
45
46 #include "bench_setup.h"
47
48 #include <optional>
49
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"
74
75 #include "bench_system.h"
76
77 namespace Nbnxm
78 {
79
80 /*! \brief Checks the kernel setup
81  *
82  * Returns an error string when the kernel is not available.
83  */
84 static std::optional<std::string> checkKernelSetup(const KernelBenchOptions& options)
85 {
86     GMX_RELEASE_ASSERT(options.nbnxmSimd < BenchMarkKernels::Count
87                                && options.nbnxmSimd != BenchMarkKernels::SimdAuto,
88                        "Need a valid kernel SIMD type");
89
90     // Check SIMD support
91     if ((options.nbnxmSimd != BenchMarkKernels::SimdNo && !GMX_SIMD)
92 #ifndef GMX_NBNXN_SIMD_4XN
93         || options.nbnxmSimd == BenchMarkKernels::Simd4XM
94 #endif
95 #ifndef GMX_NBNXN_SIMD_2XNN
96         || options.nbnxmSimd == BenchMarkKernels::Simd2XMM
97 #endif
98     )
99     {
100         return "the requested SIMD kernel was not set up at configuration time";
101     }
102
103     if (options.reportTime && (0 > gmx_cycles_calibrate(1.0)))
104     {
105         return "the -time option is not supported on this system";
106     }
107
108     return {};
109 }
110
111 //! Helper to translate between the different enumeration values.
112 static KernelType translateBenchmarkEnum(const BenchMarkKernels& kernel)
113 {
114     int kernelInt = static_cast<int>(kernel);
115     return static_cast<KernelType>(kernelInt);
116 }
117
118 /*! \brief Returns the kernel setup
119  */
120 static KernelSetup getKernelSetup(const KernelBenchOptions& options)
121 {
122     auto messageWhenInvalid = checkKernelSetup(options);
123     GMX_RELEASE_ASSERT(!messageWhenInvalid, "Need valid options");
124
125     KernelSetup kernelSetup;
126
127     // The int enum options.nbnxnSimd is set up to match KernelType + 1
128     kernelSetup.kernelType = translateBenchmarkEnum(options.nbnxmSimd);
129     // The plain-C kernel does not support analytical ewald correction
130     if (kernelSetup.kernelType == KernelType::Cpu4x4_PlainC)
131     {
132         kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
133     }
134     else
135     {
136         kernelSetup.ewaldExclusionType = options.useTabulatedEwaldCorr ? EwaldExclusionType::Table
137                                                                        : EwaldExclusionType::Analytical;
138     }
139
140     return kernelSetup;
141 }
142
143 //! Return an interaction constants struct with members used in the benchmark set appropriately
144 static interaction_const_t setupInteractionConst(const KernelBenchOptions& options)
145
146 {
147     interaction_const_t ic;
148
149     ic.vdwtype      = VanDerWaalsType::Cut;
150     ic.vdw_modifier = InteractionModifiers::PotShift;
151     ic.rvdw         = options.pairlistCutoff;
152
153     ic.eeltype = (options.coulombType == BenchMarkCoulomb::Pme ? CoulombInteractionType::Pme
154                                                                : CoulombInteractionType::RF);
155     ic.coulomb_modifier = InteractionModifiers::PotShift;
156     ic.rcoulomb         = options.pairlistCutoff;
157
158     // Reaction-field with reactionFieldPermitivity=inf
159     // TODO: Replace by calc_rffac() after refactoring that
160     ic.reactionFieldCoefficient = 0.5 * std::pow(ic.rcoulomb, -3);
161     ic.reactionFieldShift = 1 / ic.rcoulomb + ic.reactionFieldCoefficient * ic.rcoulomb * ic.rcoulomb;
162
163     if (EEL_PME_EWALD(ic.eeltype))
164     {
165         // Ewald coefficients, we ignore the potential shift
166         GMX_RELEASE_ASSERT(options.ewaldcoeff_q > 0, "Ewald coefficient should be > 0");
167         ic.ewaldcoeff_q       = options.ewaldcoeff_q;
168         ic.coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
169         init_interaction_const_tables(nullptr, &ic, 0, 0);
170     }
171
172     return ic;
173 }
174
175 //! Sets up and returns a Nbnxm object for the given benchmark options and system
176 static std::unique_ptr<nonbonded_verlet_t> setupNbnxmForBenchInstance(const KernelBenchOptions& options,
177                                                                       const gmx::BenchmarkSystem& system)
178 {
179     const auto pinPolicy  = (options.useGpu ? gmx::PinningPolicy::PinnedIfSupported
180                                            : gmx::PinningPolicy::CannotBePinned);
181     const int  numThreads = options.numThreads;
182     // Note: the options and Nbnxm combination rule enums values should match
183     const int combinationRule = static_cast<int>(options.ljCombinationRule);
184
185     auto messageWhenInvalid = checkKernelSetup(options);
186     if (messageWhenInvalid)
187     {
188         gmx_fatal(FARGS, "Requested kernel is unavailable because %s.", messageWhenInvalid->c_str());
189     }
190     Nbnxm::KernelSetup kernelSetup = getKernelSetup(options);
191
192     PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
193
194     GridSet gridSet(
195             PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy);
196
197     auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
198
199     auto pairSearch = std::make_unique<PairSearch>(
200             PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy);
201
202     auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy);
203
204     // Put everything together
205     auto nbv = std::make_unique<nonbonded_verlet_t>(
206             std::move(pairlistSets), std::move(pairSearch), std::move(atomData), kernelSetup, nullptr, nullptr);
207
208     nbnxn_atomdata_init(gmx::MDLogger(),
209                         nbv->nbat.get(),
210                         kernelSetup.kernelType,
211                         combinationRule,
212                         system.numAtomTypes,
213                         system.nonbondedParameters,
214                         1,
215                         numThreads);
216
217     t_nrnb nrnb;
218
219     GMX_RELEASE_ASSERT(!TRICLINIC(system.box), "Only rectangular unit-cells are supported here");
220     const rvec lowerCorner = { 0, 0, 0 };
221     const rvec upperCorner = { system.box[XX][XX], system.box[YY][YY], system.box[ZZ][ZZ] };
222
223     gmx::ArrayRef<const int> atomInfo;
224     if (options.useHalfLJOptimization)
225     {
226         atomInfo = system.atomInfoOxygenVdw;
227     }
228     else
229     {
230         atomInfo = system.atomInfoAllVdw;
231     }
232
233     const real atomDensity = system.coordinates.size() / det(system.box);
234
235     nbnxn_put_on_grid(nbv.get(),
236                       system.box,
237                       0,
238                       lowerCorner,
239                       upperCorner,
240                       nullptr,
241                       { 0, int(system.coordinates.size()) },
242                       atomDensity,
243                       atomInfo,
244                       system.coordinates,
245                       0,
246                       nullptr);
247
248     nbv->constructPairlist(gmx::InteractionLocality::Local, system.excls, 0, &nrnb);
249
250     nbv->setAtomProperties(system.atomTypes, system.charges, atomInfo);
251
252     return nbv;
253 }
254
255 //! Add the options instance to the list for all requested kernel SIMD types
256 static void expandSimdOptionAndPushBack(const KernelBenchOptions&        options,
257                                         std::vector<KernelBenchOptions>* optionsList)
258 {
259     if (options.nbnxmSimd == BenchMarkKernels::SimdAuto)
260     {
261         bool addedInstance = false;
262 #ifdef GMX_NBNXN_SIMD_4XN
263         optionsList->push_back(options);
264         optionsList->back().nbnxmSimd = BenchMarkKernels::Simd4XM;
265         addedInstance                 = true;
266 #endif
267 #ifdef GMX_NBNXN_SIMD_2XNN
268         optionsList->push_back(options);
269         optionsList->back().nbnxmSimd = BenchMarkKernels::Simd2XMM;
270         addedInstance                 = true;
271 #endif
272         if (!addedInstance)
273         {
274             optionsList->push_back(options);
275             optionsList->back().nbnxmSimd = BenchMarkKernels::SimdNo;
276         }
277     }
278     else
279     {
280         optionsList->push_back(options);
281     }
282 }
283
284 //! Sets up and runs the requested benchmark instance and prints the results
285 //
286 // When \p doWarmup is true runs the warmup iterations instead
287 // of the normal ones and does not print any results
288 static void setupAndRunInstance(const gmx::BenchmarkSystem& system,
289                                 const KernelBenchOptions&   options,
290                                 const bool                  doWarmup)
291 {
292     // Generate an, accurate, estimate of the number of non-zero pair interactions
293     const real atomDensity = system.coordinates.size() / det(system.box);
294     const real numPairsWithinCutoff =
295             atomDensity * 4.0 / 3.0 * M_PI * std::pow(options.pairlistCutoff, 3);
296     const real numUsefulPairs = system.coordinates.size() * 0.5 * (numPairsWithinCutoff + 1);
297
298     std::unique_ptr<nonbonded_verlet_t> nbv = setupNbnxmForBenchInstance(options, system);
299
300     // We set the interaction cut-off to the pairlist cut-off
301     interaction_const_t ic = setupInteractionConst(options);
302
303     t_nrnb nrnb = { 0 };
304
305     gmx_enerdata_t enerd(1, 0);
306
307     gmx::StepWorkload stepWork;
308     stepWork.computeForces = true;
309     if (options.computeVirialAndEnergy)
310     {
311         stepWork.computeVirial = true;
312         stepWork.computeEnergy = true;
313     }
314
315     const gmx::EnumerationArray<BenchMarkKernels, std::string> kernelNames = {
316         "auto", "no", "4xM", "2xMM"
317     };
318
319     const gmx::EnumerationArray<BenchMarkCombRule, std::string> combruleNames = { "geom.",
320                                                                                   "LB",
321                                                                                   "none" };
322
323     if (!doWarmup)
324     {
325         fprintf(stdout,
326                 "%-7s %-4s %-5s %-4s ",
327                 options.coulombType == BenchMarkCoulomb::Pme ? "Ewald" : "RF",
328                 options.useHalfLJOptimization ? "half" : "all",
329                 combruleNames[options.ljCombinationRule].c_str(),
330                 kernelNames[options.nbnxmSimd].c_str());
331         if (!options.outputFile.empty())
332         {
333             fprintf(system.csv,
334                     "\"%d\",\"%zu\",\"%g\",\"%d\",\"%d\",\"%s\",\"%s\",\"%s\",\"%s\",\"%s\",\"%"
335                     "s\",",
336 #if GMX_SIMD
337                     (options.nbnxmSimd != BenchMarkKernels::SimdNo) ? GMX_SIMD_REAL_WIDTH : 0,
338 #else
339                     0,
340 #endif
341                     system.coordinates.size(),
342                     options.pairlistCutoff,
343                     options.numThreads,
344                     options.numIterations,
345                     options.computeVirialAndEnergy ? "yes" : "no",
346                     (options.coulombType != BenchMarkCoulomb::ReactionField)
347                             ? ((options.nbnxmSimd == BenchMarkKernels::SimdNo || options.useTabulatedEwaldCorr)
348                                        ? "table"
349                                        : "analytical")
350                             : "",
351                     options.coulombType == BenchMarkCoulomb::Pme ? "Ewald" : "RF",
352                     options.useHalfLJOptimization ? "half" : "all",
353                     combruleNames[options.ljCombinationRule].c_str(),
354                     kernelNames[options.nbnxmSimd].c_str());
355         }
356     }
357
358     // Run pre-iteration to avoid cache misses
359     for (int iter = 0; iter < options.numPreIterations; iter++)
360     {
361         nbv->dispatchNonbondedKernel(
362                 gmx::InteractionLocality::Local, ic, stepWork, enbvClearFYes, system.forceRec, &enerd, &nrnb);
363     }
364
365     const int numIterations = (doWarmup ? options.numWarmupIterations : options.numIterations);
366     const PairlistSet& pairlistSet = nbv->pairlistSets().pairlistSet(gmx::InteractionLocality::Local);
367     const gmx::index numPairs = pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_ + pairlistSet.natpair_q_;
368     gmx_cycles_t cycles = gmx_cycles_read();
369     for (int iter = 0; iter < numIterations; iter++)
370     {
371         // Run the kernel without force clearing
372         nbv->dispatchNonbondedKernel(
373                 gmx::InteractionLocality::Local, ic, stepWork, enbvClearFNo, system.forceRec, &enerd, &nrnb);
374     }
375     cycles = gmx_cycles_read() - cycles;
376     if (!doWarmup)
377     {
378         if (options.reportTime)
379         {
380             const double uSec = static_cast<double>(cycles) * gmx_cycles_calibrate(1.0) * 1.e6;
381             if (options.cyclesPerPair)
382             {
383                 fprintf(stdout,
384                         "%13.2f %13.3f %10.3f %10.3f\n",
385                         uSec,
386                         uSec / options.numIterations,
387                         uSec / (options.numIterations * numPairs),
388                         uSec / (options.numIterations * numUsefulPairs));
389                 if (!options.outputFile.empty())
390                 {
391                     fprintf(system.csv,
392                             "\"%.3f\",\"%.4f\",\"%.4f\",\"%.4f\"\n",
393                             uSec,
394                             uSec / options.numIterations,
395                             uSec / (options.numIterations * numPairs),
396                             uSec / (options.numIterations * numUsefulPairs));
397                 }
398             }
399             else
400             {
401                 fprintf(stdout,
402                         "%13.2f %13.3f %10.3f %10.3f\n",
403                         uSec,
404                         uSec / options.numIterations,
405                         options.numIterations * numPairs / uSec,
406                         options.numIterations * numUsefulPairs / uSec);
407                 if (!options.outputFile.empty())
408                 {
409                     fprintf(system.csv,
410                             "\"%.3f\",\"%.4f\",\"%.4f\",\"%.4f\"\n",
411                             uSec,
412                             uSec / options.numIterations,
413                             options.numIterations * numPairs / uSec,
414                             options.numIterations * numUsefulPairs / uSec);
415                 }
416             }
417         }
418         else
419         {
420             const double dCycles = static_cast<double>(cycles);
421             if (options.cyclesPerPair)
422             {
423                 fprintf(stdout,
424                         "%10.3f %10.4f %8.4f %8.4f\n",
425                         cycles * 1e-6,
426                         dCycles / options.numIterations * 1e-6,
427                         dCycles / (options.numIterations * numPairs),
428                         dCycles / (options.numIterations * numUsefulPairs));
429             }
430             else
431             {
432                 fprintf(stdout,
433                         "%10.3f %10.4f %8.4f %8.4f\n",
434                         dCycles * 1e-6,
435                         dCycles / options.numIterations * 1e-6,
436                         options.numIterations * numPairs / dCycles,
437                         options.numIterations * numUsefulPairs / dCycles);
438             }
439         }
440     }
441 }
442
443 void bench(const int sizeFactor, const KernelBenchOptions& options)
444 {
445     // We don't want to call gmx_omp_nthreads_init(), so we init what we need
446     gmx_omp_nthreads_set(emntPairsearch, options.numThreads);
447     gmx_omp_nthreads_set(emntNonbonded, options.numThreads);
448
449     const gmx::BenchmarkSystem system(sizeFactor, options.outputFile);
450
451     real minBoxSize = norm(system.box[XX]);
452     for (int dim = YY; dim < DIM; dim++)
453     {
454         minBoxSize = std::min(minBoxSize, norm(system.box[dim]));
455     }
456     if (options.pairlistCutoff > 0.5 * minBoxSize)
457     {
458         gmx_fatal(FARGS, "The cut-off should be shorter than half the box size");
459     }
460
461     std::vector<KernelBenchOptions> optionsList;
462     if (options.doAll)
463     {
464         KernelBenchOptions                        opt = options;
465         gmx::EnumerationWrapper<BenchMarkCoulomb> coulombIter;
466         for (auto coulombType : coulombIter)
467         {
468             opt.coulombType = coulombType;
469             for (int halfLJ = 0; halfLJ <= 1; halfLJ++)
470             {
471                 opt.useHalfLJOptimization = (halfLJ == 1);
472
473                 gmx::EnumerationWrapper<BenchMarkCombRule> combRuleIter;
474                 for (auto combRule : combRuleIter)
475                 {
476                     opt.ljCombinationRule = combRule;
477
478                     expandSimdOptionAndPushBack(opt, &optionsList);
479                 }
480             }
481         }
482     }
483     else
484     {
485         expandSimdOptionAndPushBack(options, &optionsList);
486     }
487     GMX_RELEASE_ASSERT(!optionsList.empty(), "Expect at least on benchmark setup");
488
489 #if GMX_SIMD
490     if (options.nbnxmSimd != BenchMarkKernels::SimdNo)
491     {
492         fprintf(stdout, "SIMD width:           %d\n", GMX_SIMD_REAL_WIDTH);
493     }
494 #endif
495     fprintf(stdout, "System size:          %zu atoms\n", system.coordinates.size());
496     fprintf(stdout, "Cut-off radius:       %g nm\n", options.pairlistCutoff);
497     fprintf(stdout, "Number of threads:    %d\n", options.numThreads);
498     fprintf(stdout, "Number of iterations: %d\n", options.numIterations);
499     fprintf(stdout, "Compute energies:     %s\n", options.computeVirialAndEnergy ? "yes" : "no");
500     if (options.coulombType != BenchMarkCoulomb::ReactionField)
501     {
502         fprintf(stdout,
503                 "Ewald excl. corr.:    %s\n",
504                 options.nbnxmSimd == BenchMarkKernels::SimdNo || options.useTabulatedEwaldCorr
505                         ? "table"
506                         : "analytical");
507     }
508     printf("\n");
509
510     if (options.numWarmupIterations > 0)
511     {
512         setupAndRunInstance(system, optionsList[0], true);
513     }
514
515     if (options.reportTime)
516     {
517         fprintf(stdout,
518                 "Coulomb LJ   comb. SIMD       usec         usec/it.        %s\n",
519                 options.cyclesPerPair ? "usec/pair" : "pairs/usec");
520         if (!options.outputFile.empty())
521         {
522             fprintf(system.csv,
523                     "\"width\",\"atoms\",\"cut-off radius\",\"threads\",\"iter\",\"compute "
524                     "energy\",\"Ewald excl. "
525                     "corr.\",\"Coulomb\",\"LJ\",\"comb\",\"SIMD\",\"usec\",\"usec/it\",\"total "
526                     "pairs/usec\",\"useful pairs/usec\"\n");
527         }
528         fprintf(stdout,
529                 "                                                        total      useful\n");
530     }
531     else
532     {
533         fprintf(stdout,
534                 "Coulomb LJ   comb. SIMD    Mcycles  Mcycles/it.   %s\n",
535                 options.cyclesPerPair ? "cycles/pair" : "pairs/cycle");
536         if (!options.outputFile.empty())
537         {
538             fprintf(system.csv,
539                     "\"width\",\"atoms\",\"cut-off radius\",\"threads\",\"iter\",\"compute "
540                     "energy\",\"Ewald excl. "
541                     "corr.\",\"Coulomb\",\"LJ\",\"comb\",\"SIMD\",\"Mcycles\",\"Mcycles/"
542                     "it\",\"total "
543                     "total cycles/pair\",\"total cycles per useful pair\"\n");
544         }
545         fprintf(stdout, "                                                total    useful\n");
546     }
547
548     for (const auto& optionsInstance : optionsList)
549     {
550         setupAndRunInstance(system, optionsInstance, false);
551     }
552
553     if (!options.outputFile.empty())
554     {
555         fclose(system.csv);
556     }
557 }
558
559 } // namespace Nbnxm