Use ListOfLists in gmx_mtop_t and gmx_localtop_t
[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, 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 "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"
73
74 #include "bench_system.h"
75
76 namespace Nbnxm
77 {
78
79 /*! \brief Checks the kernel setup
80  *
81  * Returns an error string when the kernel is not available.
82  */
83 static gmx::compat::optional<std::string> checkKernelSetup(const KernelBenchOptions& options)
84 {
85     GMX_RELEASE_ASSERT(options.nbnxmSimd < BenchMarkKernels::Count
86                                && options.nbnxmSimd != BenchMarkKernels::SimdAuto,
87                        "Need a valid kernel SIMD type");
88
89     // Check SIMD support
90     if ((options.nbnxmSimd != BenchMarkKernels::SimdNo && !GMX_SIMD)
91 #ifndef GMX_NBNXN_SIMD_4XN
92         || options.nbnxmSimd == BenchMarkKernels::Simd4XM
93 #endif
94 #ifndef GMX_NBNXN_SIMD_2XNN
95         || options.nbnxmSimd == BenchMarkKernels::Simd2XMM
96 #endif
97     )
98     {
99         return "the requested SIMD kernel was not set up at configuration time";
100     }
101
102     return {};
103 }
104
105 //! Helper to translate between the different enumeration values.
106 static KernelType translateBenchmarkEnum(const BenchMarkKernels& kernel)
107 {
108     int kernelInt = static_cast<int>(kernel);
109     return static_cast<KernelType>(kernelInt);
110 }
111
112 /*! \brief Returns the kernel setup
113  */
114 static KernelSetup getKernelSetup(const KernelBenchOptions& options)
115 {
116     auto messageWhenInvalid = checkKernelSetup(options);
117     GMX_RELEASE_ASSERT(!messageWhenInvalid, "Need valid options");
118
119     KernelSetup kernelSetup;
120
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)
125     {
126         kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
127     }
128     else
129     {
130         kernelSetup.ewaldExclusionType = options.useTabulatedEwaldCorr ? EwaldExclusionType::Table
131                                                                        : EwaldExclusionType::Analytical;
132     }
133
134     return kernelSetup;
135 }
136
137 //! Return an interaction constants struct with members used in the benchmark set appropriately
138 static interaction_const_t setupInteractionConst(const KernelBenchOptions& options)
139
140 {
141     interaction_const_t ic;
142
143     ic.vdwtype      = evdwCUT;
144     ic.vdw_modifier = eintmodPOTSHIFT;
145     ic.rvdw         = options.pairlistCutoff;
146
147     ic.eeltype          = (options.coulombType == BenchMarkCoulomb::Pme ? eelPME : eelRF);
148     ic.coulomb_modifier = eintmodPOTSHIFT;
149     ic.rcoulomb         = options.pairlistCutoff;
150
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;
155
156     if (EEL_PME_EWALD(ic.eeltype))
157     {
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);
163     }
164
165     return ic;
166 }
167
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)
171 {
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);
177
178     auto messageWhenInvalid = checkKernelSetup(options);
179     if (messageWhenInvalid)
180     {
181         gmx_fatal(FARGS, "Requested kernel is unavailable because %s.", messageWhenInvalid->c_str());
182     }
183     Nbnxm::KernelSetup kernelSetup = getKernelSetup(options);
184
185     PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
186
187     GridSet gridSet(epbcXYZ, false, nullptr, nullptr, pairlistParams.pairlistType, false,
188                     numThreads, pinPolicy);
189
190     auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
191
192     auto pairSearch = std::make_unique<PairSearch>(
193             epbcXYZ, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy);
194
195     auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy);
196
197     // Put everything together
198     auto nbv = std::make_unique<nonbonded_verlet_t>(std::move(pairlistSets), std::move(pairSearch),
199                                                     std::move(atomData), kernelSetup, nullptr, nullptr);
200
201     nbnxn_atomdata_init(gmx::MDLogger(), nbv->nbat.get(), kernelSetup.kernelType, combinationRule,
202                         system.numAtomTypes, system.nonbondedParameters, 1, numThreads);
203
204     t_nrnb nrnb;
205
206     GMX_RELEASE_ASSERT(!TRICLINIC(system.box), "Only rectangular unit-cells are supported here");
207     const rvec lowerCorner = { 0, 0, 0 };
208     const rvec upperCorner = { system.box[XX][XX], system.box[YY][YY], system.box[ZZ][ZZ] };
209
210     gmx::ArrayRef<const int> atomInfo;
211     if (options.useHalfLJOptimization)
212     {
213         atomInfo = system.atomInfoOxygenVdw;
214     }
215     else
216     {
217         atomInfo = system.atomInfoAllVdw;
218     }
219
220     const real atomDensity = system.coordinates.size() / det(system.box);
221
222     nbnxn_put_on_grid(nbv.get(), system.box, 0, lowerCorner, upperCorner, nullptr,
223                       { 0, int(system.coordinates.size()) }, atomDensity, atomInfo,
224                       system.coordinates, 0, nullptr);
225
226     nbv->constructPairlist(gmx::InteractionLocality::Local, system.excls, 0, &nrnb);
227
228     t_mdatoms mdatoms;
229     // We only use (read) the atom type and charge from mdatoms
230     mdatoms.typeA   = const_cast<int*>(system.atomTypes.data());
231     mdatoms.chargeA = const_cast<real*>(system.charges.data());
232     nbv->setAtomProperties(mdatoms, atomInfo);
233
234     return nbv;
235 }
236
237 //! Add the options instance to the list for all requested kernel SIMD types
238 static void expandSimdOptionAndPushBack(const KernelBenchOptions&        options,
239                                         std::vector<KernelBenchOptions>* optionsList)
240 {
241     if (options.nbnxmSimd == BenchMarkKernels::SimdAuto)
242     {
243         bool addedInstance = false;
244 #ifdef GMX_NBNXN_SIMD_4XN
245         optionsList->push_back(options);
246         optionsList->back().nbnxmSimd = BenchMarkKernels::Simd4XM;
247         addedInstance                 = true;
248 #endif
249 #ifdef GMX_NBNXN_SIMD_2XNN
250         optionsList->push_back(options);
251         optionsList->back().nbnxmSimd = BenchMarkKernels::Simd2XMM;
252         addedInstance                 = true;
253 #endif
254         if (!addedInstance)
255         {
256             optionsList->push_back(options);
257             optionsList->back().nbnxmSimd = BenchMarkKernels::SimdNo;
258         }
259     }
260     else
261     {
262         optionsList->push_back(options);
263     }
264 }
265
266 //! Sets up and runs the requested benchmark instance and prints the results
267 //
268 // When \p doWarmup is true runs the warmup iterations instead
269 // of the normal ones and does not print any results
270 static void setupAndRunInstance(const gmx::BenchmarkSystem& system,
271                                 const KernelBenchOptions&   options,
272                                 const bool                  doWarmup)
273 {
274     // Generate an, accurate, estimate of the number of non-zero pair interactions
275     const real atomDensity = system.coordinates.size() / det(system.box);
276     const real numPairsWithinCutoff =
277             atomDensity * 4.0 / 3.0 * M_PI * std::pow(options.pairlistCutoff, 3);
278     const real numUsefulPairs = system.coordinates.size() * 0.5 * (numPairsWithinCutoff + 1);
279
280     std::unique_ptr<nonbonded_verlet_t> nbv = setupNbnxmForBenchInstance(options, system);
281
282     // We set the interaction cut-off to the pairlist cut-off
283     interaction_const_t ic = setupInteractionConst(options);
284
285     t_nrnb nrnb = { 0 };
286
287     gmx_enerdata_t enerd(1, 0);
288
289     gmx::StepWorkload stepWork;
290     stepWork.computeForces = true;
291     if (options.computeVirialAndEnergy)
292     {
293         stepWork.computeVirial = true;
294         stepWork.computeEnergy = true;
295     }
296
297     const gmx::EnumerationArray<BenchMarkKernels, std::string> kernelNames = { "auto", "no", "4xM",
298                                                                                "2xMM" };
299
300     const gmx::EnumerationArray<BenchMarkCombRule, std::string> combruleNames = { "geom.", "LB",
301                                                                                   "none" };
302
303     if (!doWarmup)
304     {
305         fprintf(stdout, "%-7s %-4s %-5s %-4s ",
306                 options.coulombType == BenchMarkCoulomb::Pme ? "Ewald" : "RF",
307                 options.useHalfLJOptimization ? "half" : "all",
308                 combruleNames[options.ljCombinationRule].c_str(), kernelNames[options.nbnxmSimd].c_str());
309     }
310
311     // Run pre-iteration to avoid cache misses
312     for (int iter = 0; iter < options.numPreIterations; iter++)
313     {
314         nbv->dispatchNonbondedKernel(gmx::InteractionLocality::Local, ic, stepWork, enbvClearFYes,
315                                      system.forceRec, &enerd, &nrnb);
316     }
317
318     const int numIterations = (doWarmup ? options.numWarmupIterations : options.numIterations);
319     const PairlistSet& pairlistSet = nbv->pairlistSets().pairlistSet(gmx::InteractionLocality::Local);
320     const gmx::index numPairs = pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_ + pairlistSet.natpair_q_;
321     gmx_cycles_t cycles = gmx_cycles_read();
322     for (int iter = 0; iter < numIterations; iter++)
323     {
324         // Run the kernel without force clearing
325         nbv->dispatchNonbondedKernel(gmx::InteractionLocality::Local, ic, stepWork, enbvClearFNo,
326                                      system.forceRec, &enerd, &nrnb);
327     }
328     cycles = gmx_cycles_read() - cycles;
329     if (!doWarmup)
330     {
331         const double dCycles = static_cast<double>(cycles);
332         if (options.cyclesPerPair)
333         {
334             fprintf(stdout, "%10.3f %10.4f %8.4f %8.4f\n", cycles * 1e-6,
335                     dCycles / options.numIterations * 1e-6, dCycles / (options.numIterations * numPairs),
336                     dCycles / (options.numIterations * numUsefulPairs));
337         }
338         else
339         {
340             fprintf(stdout, "%10.3f %10.4f %8.4f %8.4f\n", dCycles * 1e-6,
341                     dCycles / options.numIterations * 1e-6, options.numIterations * numPairs / dCycles,
342                     options.numIterations * numUsefulPairs / dCycles);
343         }
344     }
345 }
346
347 void bench(const int sizeFactor, const KernelBenchOptions& options)
348 {
349     // We don't want to call gmx_omp_nthreads_init(), so we init what we need
350     gmx_omp_nthreads_set(emntPairsearch, options.numThreads);
351     gmx_omp_nthreads_set(emntNonbonded, options.numThreads);
352
353     const gmx::BenchmarkSystem system(sizeFactor);
354
355     real minBoxSize = norm(system.box[XX]);
356     for (int dim = YY; dim < DIM; dim++)
357     {
358         minBoxSize = std::min(minBoxSize, norm(system.box[dim]));
359     }
360     if (options.pairlistCutoff > 0.5 * minBoxSize)
361     {
362         gmx_fatal(FARGS, "The cut-off should be shorter than half the box size");
363     }
364
365     std::vector<KernelBenchOptions> optionsList;
366     if (options.doAll)
367     {
368         KernelBenchOptions                        opt = options;
369         gmx::EnumerationWrapper<BenchMarkCoulomb> coulombIter;
370         for (auto coulombType : coulombIter)
371         {
372             opt.coulombType = coulombType;
373             for (int halfLJ = 0; halfLJ <= 1; halfLJ++)
374             {
375                 opt.useHalfLJOptimization = (halfLJ == 1);
376
377                 gmx::EnumerationWrapper<BenchMarkCombRule> combRuleIter;
378                 for (auto combRule : combRuleIter)
379                 {
380                     if (combRule == BenchMarkCombRule::RuleNone)
381                     {
382                         continue;
383                     }
384                     opt.ljCombinationRule = combRule;
385
386                     expandSimdOptionAndPushBack(opt, &optionsList);
387                 }
388             }
389         }
390     }
391     else
392     {
393         expandSimdOptionAndPushBack(options, &optionsList);
394     }
395     GMX_RELEASE_ASSERT(!optionsList.empty(), "Expect at least on benchmark setup");
396
397 #if GMX_SIMD
398     if (options.nbnxmSimd != BenchMarkKernels::SimdNo)
399     {
400         fprintf(stdout, "SIMD width:           %d\n", GMX_SIMD_REAL_WIDTH);
401     }
402 #endif
403     fprintf(stdout, "System size:          %zu atoms\n", system.coordinates.size());
404     fprintf(stdout, "Cut-off radius:       %g nm\n", options.pairlistCutoff);
405     fprintf(stdout, "Number of threads:    %d\n", options.numThreads);
406     fprintf(stdout, "Number of iterations: %d\n", options.numIterations);
407     fprintf(stdout, "Compute energies:     %s\n", options.computeVirialAndEnergy ? "yes" : "no");
408     if (options.coulombType != BenchMarkCoulomb::ReactionField)
409     {
410         fprintf(stdout, "Ewald excl. corr.:    %s\n",
411                 options.nbnxmSimd == BenchMarkKernels::SimdNo || options.useTabulatedEwaldCorr
412                         ? "table"
413                         : "analytical");
414     }
415     printf("\n");
416
417     if (options.numWarmupIterations > 0)
418     {
419         setupAndRunInstance(system, optionsList[0], true);
420     }
421
422     fprintf(stdout, "Coulomb LJ   comb. SIMD    Mcycles  Mcycles/it.   %s\n",
423             options.cyclesPerPair ? "cycles/pair" : "pairs/cycle");
424     fprintf(stdout, "                                                total    useful\n");
425
426     for (const auto& optionsInstance : optionsList)
427     {
428         setupAndRunInstance(system, optionsInstance, false);
429     }
430 }
431
432 } // namespace Nbnxm