A check for perturbed listed pairs beyond rlist
[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, 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      = evdwCUT;
150     ic.vdw_modifier = eintmodPOTSHIFT;
151     ic.rvdw         = options.pairlistCutoff;
152
153     ic.eeltype          = (options.coulombType == BenchMarkCoulomb::Pme ? eelPME : eelRF);
154     ic.coulomb_modifier = eintmodPOTSHIFT;
155     ic.rcoulomb         = options.pairlistCutoff;
156
157     // Reaction-field with epsilon_rf=inf
158     // TODO: Replace by calc_rffac() after refactoring that
159     ic.k_rf = 0.5 * std::pow(ic.rcoulomb, -3);
160     ic.c_rf = 1 / ic.rcoulomb + ic.k_rf * ic.rcoulomb * ic.rcoulomb;
161
162     if (EEL_PME_EWALD(ic.eeltype))
163     {
164         // Ewald coefficients, we ignore the potential shift
165         GMX_RELEASE_ASSERT(options.ewaldcoeff_q > 0, "Ewald coefficient should be > 0");
166         ic.ewaldcoeff_q       = options.ewaldcoeff_q;
167         ic.coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
168         init_interaction_const_tables(nullptr, &ic, 0, 0);
169     }
170
171     return ic;
172 }
173
174 //! Sets up and returns a Nbnxm object for the given benchmark options and system
175 static std::unique_ptr<nonbonded_verlet_t> setupNbnxmForBenchInstance(const KernelBenchOptions& options,
176                                                                       const gmx::BenchmarkSystem& system)
177 {
178     const auto pinPolicy  = (options.useGpu ? gmx::PinningPolicy::PinnedIfSupported
179                                            : gmx::PinningPolicy::CannotBePinned);
180     const int  numThreads = options.numThreads;
181     // Note: the options and Nbnxm combination rule enums values should match
182     const int combinationRule = static_cast<int>(options.ljCombinationRule);
183
184     auto messageWhenInvalid = checkKernelSetup(options);
185     if (messageWhenInvalid)
186     {
187         gmx_fatal(FARGS, "Requested kernel is unavailable because %s.", messageWhenInvalid->c_str());
188     }
189     Nbnxm::KernelSetup kernelSetup = getKernelSetup(options);
190
191     PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
192
193     GridSet gridSet(PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false,
194                     numThreads, pinPolicy);
195
196     auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
197
198     auto pairSearch =
199             std::make_unique<PairSearch>(PbcType::Xyz, false, nullptr, nullptr,
200                                          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>(std::move(pairlistSets), std::move(pairSearch),
206                                                     std::move(atomData), kernelSetup, nullptr, nullptr);
207
208     nbnxn_atomdata_init(gmx::MDLogger(), nbv->nbat.get(), kernelSetup.kernelType, combinationRule,
209                         system.numAtomTypes, system.nonbondedParameters, 1, numThreads);
210
211     t_nrnb nrnb;
212
213     GMX_RELEASE_ASSERT(!TRICLINIC(system.box), "Only rectangular unit-cells are supported here");
214     const rvec lowerCorner = { 0, 0, 0 };
215     const rvec upperCorner = { system.box[XX][XX], system.box[YY][YY], system.box[ZZ][ZZ] };
216
217     gmx::ArrayRef<const int> atomInfo;
218     if (options.useHalfLJOptimization)
219     {
220         atomInfo = system.atomInfoOxygenVdw;
221     }
222     else
223     {
224         atomInfo = system.atomInfoAllVdw;
225     }
226
227     const real atomDensity = system.coordinates.size() / det(system.box);
228
229     nbnxn_put_on_grid(nbv.get(), system.box, 0, lowerCorner, upperCorner, nullptr,
230                       { 0, int(system.coordinates.size()) }, atomDensity, atomInfo,
231                       system.coordinates, 0, nullptr);
232
233     nbv->constructPairlist(gmx::InteractionLocality::Local, system.excls, 0, &nrnb);
234
235     nbv->setAtomProperties(system.atomTypes, system.charges, atomInfo);
236
237     return nbv;
238 }
239
240 //! Add the options instance to the list for all requested kernel SIMD types
241 static void expandSimdOptionAndPushBack(const KernelBenchOptions&        options,
242                                         std::vector<KernelBenchOptions>* optionsList)
243 {
244     if (options.nbnxmSimd == BenchMarkKernels::SimdAuto)
245     {
246         bool addedInstance = false;
247 #ifdef GMX_NBNXN_SIMD_4XN
248         optionsList->push_back(options);
249         optionsList->back().nbnxmSimd = BenchMarkKernels::Simd4XM;
250         addedInstance                 = true;
251 #endif
252 #ifdef GMX_NBNXN_SIMD_2XNN
253         optionsList->push_back(options);
254         optionsList->back().nbnxmSimd = BenchMarkKernels::Simd2XMM;
255         addedInstance                 = true;
256 #endif
257         if (!addedInstance)
258         {
259             optionsList->push_back(options);
260             optionsList->back().nbnxmSimd = BenchMarkKernels::SimdNo;
261         }
262     }
263     else
264     {
265         optionsList->push_back(options);
266     }
267 }
268
269 //! Sets up and runs the requested benchmark instance and prints the results
270 //
271 // When \p doWarmup is true runs the warmup iterations instead
272 // of the normal ones and does not print any results
273 static void setupAndRunInstance(const gmx::BenchmarkSystem& system,
274                                 const KernelBenchOptions&   options,
275                                 const bool                  doWarmup)
276 {
277     // Generate an, accurate, estimate of the number of non-zero pair interactions
278     const real atomDensity = system.coordinates.size() / det(system.box);
279     const real numPairsWithinCutoff =
280             atomDensity * 4.0 / 3.0 * M_PI * std::pow(options.pairlistCutoff, 3);
281     const real numUsefulPairs = system.coordinates.size() * 0.5 * (numPairsWithinCutoff + 1);
282
283     std::unique_ptr<nonbonded_verlet_t> nbv = setupNbnxmForBenchInstance(options, system);
284
285     // We set the interaction cut-off to the pairlist cut-off
286     interaction_const_t ic = setupInteractionConst(options);
287
288     t_nrnb nrnb = { 0 };
289
290     gmx_enerdata_t enerd(1, 0);
291
292     gmx::StepWorkload stepWork;
293     stepWork.computeForces = true;
294     if (options.computeVirialAndEnergy)
295     {
296         stepWork.computeVirial = true;
297         stepWork.computeEnergy = true;
298     }
299
300     const gmx::EnumerationArray<BenchMarkKernels, std::string> kernelNames = { "auto", "no", "4xM",
301                                                                                "2xMM" };
302
303     const gmx::EnumerationArray<BenchMarkCombRule, std::string> combruleNames = { "geom.", "LB",
304                                                                                   "none" };
305
306     if (!doWarmup)
307     {
308         fprintf(stdout, "%-7s %-4s %-5s %-4s ",
309                 options.coulombType == BenchMarkCoulomb::Pme ? "Ewald" : "RF",
310                 options.useHalfLJOptimization ? "half" : "all",
311                 combruleNames[options.ljCombinationRule].c_str(), kernelNames[options.nbnxmSimd].c_str());
312         if (!options.outputFile.empty())
313         {
314             fprintf(system.csv,
315                     "\"%d\",\"%zu\",\"%g\",\"%d\",\"%d\",\"%s\",\"%s\",\"%s\",\"%s\",\"%s\",\"%"
316                     "s\",",
317 #if GMX_SIMD
318                     (options.nbnxmSimd != BenchMarkKernels::SimdNo) ? GMX_SIMD_REAL_WIDTH : 0,
319 #else
320                     0,
321 #endif
322                     system.coordinates.size(), options.pairlistCutoff, options.numThreads,
323                     options.numIterations, options.computeVirialAndEnergy ? "yes" : "no",
324                     (options.coulombType != BenchMarkCoulomb::ReactionField)
325                             ? ((options.nbnxmSimd == BenchMarkKernels::SimdNo || options.useTabulatedEwaldCorr)
326                                        ? "table"
327                                        : "analytical")
328                             : "",
329                     options.coulombType == BenchMarkCoulomb::Pme ? "Ewald" : "RF",
330                     options.useHalfLJOptimization ? "half" : "all",
331                     combruleNames[options.ljCombinationRule].c_str(),
332                     kernelNames[options.nbnxmSimd].c_str());
333         }
334     }
335
336     // Run pre-iteration to avoid cache misses
337     for (int iter = 0; iter < options.numPreIterations; iter++)
338     {
339         nbv->dispatchNonbondedKernel(gmx::InteractionLocality::Local, ic, stepWork, enbvClearFYes,
340                                      system.forceRec, &enerd, &nrnb);
341     }
342
343     const int numIterations = (doWarmup ? options.numWarmupIterations : options.numIterations);
344     const PairlistSet& pairlistSet = nbv->pairlistSets().pairlistSet(gmx::InteractionLocality::Local);
345     const gmx::index numPairs = pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_ + pairlistSet.natpair_q_;
346     gmx_cycles_t cycles = gmx_cycles_read();
347     for (int iter = 0; iter < numIterations; iter++)
348     {
349         // Run the kernel without force clearing
350         nbv->dispatchNonbondedKernel(gmx::InteractionLocality::Local, ic, stepWork, enbvClearFNo,
351                                      system.forceRec, &enerd, &nrnb);
352     }
353     cycles = gmx_cycles_read() - cycles;
354     if (!doWarmup)
355     {
356         if (options.reportTime)
357         {
358             const double uSec = static_cast<double>(cycles) * gmx_cycles_calibrate(1.0) * 1.e6;
359             if (options.cyclesPerPair)
360             {
361                 fprintf(stdout, "%13.2f %13.3f %10.3f %10.3f\n", uSec, uSec / options.numIterations,
362                         uSec / (options.numIterations * numPairs),
363                         uSec / (options.numIterations * numUsefulPairs));
364                 if (!options.outputFile.empty())
365                 {
366                     fprintf(system.csv, "\"%.3f\",\"%.4f\",\"%.4f\",\"%.4f\"\n", uSec,
367                             uSec / options.numIterations, uSec / (options.numIterations * numPairs),
368                             uSec / (options.numIterations * numUsefulPairs));
369                 }
370             }
371             else
372             {
373                 fprintf(stdout, "%13.2f %13.3f %10.3f %10.3f\n", uSec, uSec / options.numIterations,
374                         options.numIterations * numPairs / uSec,
375                         options.numIterations * numUsefulPairs / uSec);
376                 if (!options.outputFile.empty())
377                 {
378                     fprintf(system.csv, "\"%.3f\",\"%.4f\",\"%.4f\",\"%.4f\"\n", uSec,
379                             uSec / options.numIterations, options.numIterations * numPairs / uSec,
380                             options.numIterations * numUsefulPairs / uSec);
381                 }
382             }
383         }
384         else
385         {
386             const double dCycles = static_cast<double>(cycles);
387             if (options.cyclesPerPair)
388             {
389                 fprintf(stdout, "%10.3f %10.4f %8.4f %8.4f\n", cycles * 1e-6,
390                         dCycles / options.numIterations * 1e-6,
391                         dCycles / (options.numIterations * numPairs),
392                         dCycles / (options.numIterations * numUsefulPairs));
393             }
394             else
395             {
396                 fprintf(stdout, "%10.3f %10.4f %8.4f %8.4f\n", dCycles * 1e-6,
397                         dCycles / options.numIterations * 1e-6, options.numIterations * numPairs / dCycles,
398                         options.numIterations * numUsefulPairs / dCycles);
399             }
400         }
401     }
402 }
403
404 void bench(const int sizeFactor, const KernelBenchOptions& options)
405 {
406     // We don't want to call gmx_omp_nthreads_init(), so we init what we need
407     gmx_omp_nthreads_set(emntPairsearch, options.numThreads);
408     gmx_omp_nthreads_set(emntNonbonded, options.numThreads);
409
410     const gmx::BenchmarkSystem system(sizeFactor, options.outputFile);
411
412     real minBoxSize = norm(system.box[XX]);
413     for (int dim = YY; dim < DIM; dim++)
414     {
415         minBoxSize = std::min(minBoxSize, norm(system.box[dim]));
416     }
417     if (options.pairlistCutoff > 0.5 * minBoxSize)
418     {
419         gmx_fatal(FARGS, "The cut-off should be shorter than half the box size");
420     }
421
422     std::vector<KernelBenchOptions> optionsList;
423     if (options.doAll)
424     {
425         KernelBenchOptions                        opt = options;
426         gmx::EnumerationWrapper<BenchMarkCoulomb> coulombIter;
427         for (auto coulombType : coulombIter)
428         {
429             opt.coulombType = coulombType;
430             for (int halfLJ = 0; halfLJ <= 1; halfLJ++)
431             {
432                 opt.useHalfLJOptimization = (halfLJ == 1);
433
434                 gmx::EnumerationWrapper<BenchMarkCombRule> combRuleIter;
435                 for (auto combRule : combRuleIter)
436                 {
437                     opt.ljCombinationRule = combRule;
438
439                     expandSimdOptionAndPushBack(opt, &optionsList);
440                 }
441             }
442         }
443     }
444     else
445     {
446         expandSimdOptionAndPushBack(options, &optionsList);
447     }
448     GMX_RELEASE_ASSERT(!optionsList.empty(), "Expect at least on benchmark setup");
449
450 #if GMX_SIMD
451     if (options.nbnxmSimd != BenchMarkKernels::SimdNo)
452     {
453         fprintf(stdout, "SIMD width:           %d\n", GMX_SIMD_REAL_WIDTH);
454     }
455 #endif
456     fprintf(stdout, "System size:          %zu atoms\n", system.coordinates.size());
457     fprintf(stdout, "Cut-off radius:       %g nm\n", options.pairlistCutoff);
458     fprintf(stdout, "Number of threads:    %d\n", options.numThreads);
459     fprintf(stdout, "Number of iterations: %d\n", options.numIterations);
460     fprintf(stdout, "Compute energies:     %s\n", options.computeVirialAndEnergy ? "yes" : "no");
461     if (options.coulombType != BenchMarkCoulomb::ReactionField)
462     {
463         fprintf(stdout, "Ewald excl. corr.:    %s\n",
464                 options.nbnxmSimd == BenchMarkKernels::SimdNo || options.useTabulatedEwaldCorr
465                         ? "table"
466                         : "analytical");
467     }
468     printf("\n");
469
470     if (options.numWarmupIterations > 0)
471     {
472         setupAndRunInstance(system, optionsList[0], true);
473     }
474
475     if (options.reportTime)
476     {
477         fprintf(stdout, "Coulomb LJ   comb. SIMD       usec         usec/it.        %s\n",
478                 options.cyclesPerPair ? "usec/pair" : "pairs/usec");
479         if (!options.outputFile.empty())
480         {
481             fprintf(system.csv,
482                     "\"width\",\"atoms\",\"cut-off radius\",\"threads\",\"iter\",\"compute "
483                     "energy\",\"Ewald excl. "
484                     "corr.\",\"Coulomb\",\"LJ\",\"comb\",\"SIMD\",\"usec\",\"usec/it\",\"total "
485                     "pairs/usec\",\"useful pairs/usec\"\n");
486         }
487         fprintf(stdout,
488                 "                                                        total      useful\n");
489     }
490     else
491     {
492         fprintf(stdout, "Coulomb LJ   comb. SIMD    Mcycles  Mcycles/it.   %s\n",
493                 options.cyclesPerPair ? "cycles/pair" : "pairs/cycle");
494         if (!options.outputFile.empty())
495         {
496             fprintf(system.csv,
497                     "\"width\",\"atoms\",\"cut-off radius\",\"threads\",\"iter\",\"compute "
498                     "energy\",\"Ewald excl. "
499                     "corr.\",\"Coulomb\",\"LJ\",\"comb\",\"SIMD\",\"Mcycles\",\"Mcycles/"
500                     "it\",\"total "
501                     "total cycles/pair\",\"total cycles per useful pair\"\n");
502         }
503         fprintf(stdout, "                                                total    useful\n");
504     }
505
506     for (const auto& optionsInstance : optionsList)
507     {
508         setupAndRunInstance(system, optionsInstance, false);
509     }
510
511     if (!options.outputFile.empty())
512     {
513         fclose(system.csv);
514     }
515 }
516
517 } // namespace Nbnxm