b279596730bb9f142cefbc613ceff9937736865c
[alexxy/gromacs.git] / src / gromacs / nbnxm / nbnxm_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 /*! \internal \file
36  * \brief Common functions for the different NBNXN GPU implementations.
37  *
38  * \author Berk Hess <hess@kth.se>
39  *
40  * \ingroup module_nbnxm
41  */
42
43 #include "gmxpre.h"
44
45 #include <memory>
46
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/hardware/hw_info.h"
50 #include "gromacs/mdlib/gmx_omp_nthreads.h"
51 #include "gromacs/mdtypes/commrec.h"
52 #include "gromacs/mdtypes/enerdata.h"
53 #include "gromacs/mdtypes/forcerec.h"
54 #include "gromacs/mdtypes/inputrec.h"
55 #include "gromacs/mdtypes/interaction_const.h"
56 #include "gromacs/nbnxm/atomdata.h"
57 #include "gromacs/nbnxm/gpu_data_mgmt.h"
58 #include "gromacs/nbnxm/nbnxm.h"
59 #include "gromacs/nbnxm/pairlist_tuning.h"
60 #include "gromacs/simd/simd.h"
61 #include "gromacs/topology/mtop_util.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/logger.h"
64
65 #include "grid.h"
66 #include "nbnxm_geometry.h"
67 #include "nbnxm_simd.h"
68 #include "pairlist.h"
69 #include "pairlistset.h"
70 #include "pairlistsets.h"
71 #include "pairsearch.h"
72
73 namespace Nbnxm
74 {
75
76 /*! \brief Resources that can be used to execute non-bonded kernels on */
77 enum class NonbondedResource : int
78 {
79     Cpu,
80     Gpu,
81     EmulateGpu
82 };
83
84 /*! \brief Returns whether CPU SIMD support exists for the given inputrec
85  *
86  * If the return value is FALSE and fplog/cr != NULL, prints a fallback
87  * message to fplog/stderr.
88  */
89 static bool nbnxn_simd_supported(const gmx::MDLogger& mdlog, const t_inputrec& inputrec)
90 {
91     if (inputrec.vdwtype == VanDerWaalsType::Pme && inputrec.ljpme_combination_rule == LongRangeVdW::LB)
92     {
93         /* LJ PME with LB combination rule does 7 mesh operations.
94          * This so slow that we don't compile SIMD non-bonded kernels
95          * for that. */
96         GMX_LOG(mdlog.warning)
97                 .asParagraph()
98                 .appendText(
99                         "LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling "
100                         "back to plain C kernels");
101         return FALSE;
102     }
103
104     return TRUE;
105 }
106
107 /*! \brief Returns the most suitable CPU kernel type and Ewald handling */
108 static KernelSetup pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused& inputrec,
109                                          const gmx_hw_info_t gmx_unused& hardwareInfo)
110 {
111     KernelSetup kernelSetup;
112
113     if (!GMX_SIMD)
114     {
115         kernelSetup.kernelType         = KernelType::Cpu4x4_PlainC;
116         kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
117     }
118     else
119     {
120 #ifdef GMX_NBNXN_SIMD_4XN
121         kernelSetup.kernelType = KernelType::Cpu4xN_Simd_4xN;
122 #endif
123 #ifdef GMX_NBNXN_SIMD_2XNN
124         kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
125 #endif
126
127 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
128         /* We need to choose if we want 2x(N+N) or 4xN kernels.
129          * This is based on the SIMD acceleration choice and CPU information
130          * detected at runtime.
131          *
132          * 4xN calculates more (zero) interactions, but has less pair-search
133          * work and much better kernel instruction scheduling.
134          *
135          * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
136          * which doesn't have FMA, both the analytical and tabulated Ewald
137          * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
138          * 2x(4+4) because it results in significantly fewer pairs.
139          * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
140          * 10% with HT, 50% without HT. As we currently don't detect the actual
141          * use of HT, use 4x8 to avoid a potential performance hit.
142          * On Intel Haswell 4x8 is always faster.
143          */
144         kernelSetup.kernelType = KernelType::Cpu4xN_Simd_4xN;
145
146         if (!GMX_SIMD_HAVE_FMA && (EEL_PME_EWALD(inputrec.coulombtype) || EVDW_PME(inputrec.vdwtype)))
147         {
148             /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
149              * There are enough instructions to make 2x(4+4) efficient.
150              */
151             kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
152         }
153
154         if (hardwareInfo.haveAmdZen1Cpu)
155         {
156             /* One 256-bit FMA per cycle makes 2xNN faster */
157             kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
158         }
159 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
160
161
162         if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
163         {
164 #ifdef GMX_NBNXN_SIMD_4XN
165             kernelSetup.kernelType = KernelType::Cpu4xN_Simd_4xN;
166 #else
167             gmx_fatal(FARGS,
168                       "SIMD 4xN kernels requested, but GROMACS has been compiled without support "
169                       "for these kernels");
170 #endif
171         }
172         if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
173         {
174 #ifdef GMX_NBNXN_SIMD_2XNN
175             kernelSetup.kernelType = KernelType::Cpu4xN_Simd_2xNN;
176 #else
177             gmx_fatal(FARGS,
178                       "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without "
179                       "support for these kernels");
180 #endif
181         }
182
183         /* Analytical Ewald exclusion correction is only an option in
184          * the SIMD kernel.
185          * Since table lookup's don't parallelize with SIMD, analytical
186          * will probably always be faster for a SIMD width of 8 or more.
187          * With FMA analytical is sometimes faster for a width if 4 as well.
188          * In single precision, this is faster on Bulldozer.
189          * On AMD Zen, tabulated Ewald kernels are faster on all 4 combinations
190          * of single or double precision and 128 or 256-bit AVX2.
191          */
192         MSVC_DIAGNOSTIC_IGNORE(6285) // Always zero because compile time constant
193         if (
194 #if GMX_SIMD
195                 (GMX_SIMD_REAL_WIDTH >= 8 || (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)) &&
196 #endif
197                 !hardwareInfo.haveAmdZen1Cpu)
198         {
199             kernelSetup.ewaldExclusionType = EwaldExclusionType::Analytical;
200         }
201         MSVC_DIAGNOSTIC_RESET
202         else { kernelSetup.ewaldExclusionType = EwaldExclusionType::Table; }
203         if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
204         {
205             kernelSetup.ewaldExclusionType = EwaldExclusionType::Table;
206         }
207         if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
208         {
209             kernelSetup.ewaldExclusionType = EwaldExclusionType::Analytical;
210         }
211     }
212
213     return kernelSetup;
214 }
215
216 const char* lookup_kernel_name(const KernelType kernelType)
217 {
218     switch (kernelType)
219     {
220         case KernelType::NotSet: return "not set";
221         case KernelType::Cpu4x4_PlainC: return "plain C";
222         case KernelType::Cpu4xN_Simd_4xN:
223         case KernelType::Cpu4xN_Simd_2xNN:
224 #if GMX_SIMD
225             return "SIMD";
226 #else  // GMX_SIMD
227             return "not available";
228 #endif // GMX_SIMD
229         case KernelType::Gpu8x8x8: return "GPU";
230         case KernelType::Cpu8x8x8_PlainC: return "plain C";
231
232         default: gmx_fatal(FARGS, "Illegal kernel type selected");
233     }
234 };
235
236 /*! \brief Returns the most suitable kernel type and Ewald handling */
237 static KernelSetup pick_nbnxn_kernel(const gmx::MDLogger&     mdlog,
238                                      gmx_bool                 use_simd_kernels,
239                                      const gmx_hw_info_t&     hardwareInfo,
240                                      const NonbondedResource& nonbondedResource,
241                                      const t_inputrec&        inputrec)
242 {
243     KernelSetup kernelSetup;
244
245     if (nonbondedResource == NonbondedResource::EmulateGpu)
246     {
247         kernelSetup.kernelType         = KernelType::Cpu8x8x8_PlainC;
248         kernelSetup.ewaldExclusionType = EwaldExclusionType::DecidedByGpuModule;
249
250         GMX_LOG(mdlog.warning).asParagraph().appendText("Emulating a GPU run on the CPU (slow)");
251     }
252     else if (nonbondedResource == NonbondedResource::Gpu)
253     {
254         kernelSetup.kernelType         = KernelType::Gpu8x8x8;
255         kernelSetup.ewaldExclusionType = EwaldExclusionType::DecidedByGpuModule;
256     }
257     else
258     {
259         if (use_simd_kernels && nbnxn_simd_supported(mdlog, inputrec))
260         {
261             kernelSetup = pick_nbnxn_kernel_cpu(inputrec, hardwareInfo);
262         }
263         else
264         {
265             kernelSetup.kernelType         = KernelType::Cpu4x4_PlainC;
266             kernelSetup.ewaldExclusionType = EwaldExclusionType::Analytical;
267         }
268     }
269
270     GMX_LOG(mdlog.info)
271             .asParagraph()
272             .appendTextFormatted("Using %s %dx%d nonbonded short-range kernels",
273                                  lookup_kernel_name(kernelSetup.kernelType),
274                                  IClusterSizePerKernelType[kernelSetup.kernelType],
275                                  JClusterSizePerKernelType[kernelSetup.kernelType]);
276
277     if (KernelType::Cpu4x4_PlainC == kernelSetup.kernelType
278         || KernelType::Cpu8x8x8_PlainC == kernelSetup.kernelType)
279     {
280         GMX_LOG(mdlog.warning)
281                 .asParagraph()
282                 .appendTextFormatted(
283                         "WARNING: Using the slow %s kernels. This should\n"
284                         "not happen during routine usage on supported platforms.",
285                         lookup_kernel_name(kernelSetup.kernelType));
286     }
287
288     GMX_RELEASE_ASSERT(kernelSetup.kernelType != KernelType::NotSet
289                                && kernelSetup.ewaldExclusionType != EwaldExclusionType::NotSet,
290                        "All kernel setup parameters should be set here");
291
292     return kernelSetup;
293 }
294
295 } // namespace Nbnxm
296
297 PairlistSets::PairlistSets(const PairlistParams& pairlistParams,
298                            const bool            haveMultipleDomains,
299                            const int             minimumIlistCountForGpuBalancing) :
300     params_(pairlistParams), minimumIlistCountForGpuBalancing_(minimumIlistCountForGpuBalancing)
301 {
302     localSet_ = std::make_unique<PairlistSet>(params_);
303
304     if (haveMultipleDomains)
305     {
306         nonlocalSet_ = std::make_unique<PairlistSet>(params_);
307     }
308 }
309
310 namespace Nbnxm
311 {
312
313 /*! \brief Gets and returns the minimum i-list count for balacing based on the GPU used or env.var. when set */
314 static int getMinimumIlistCountForGpuBalancing(NbnxmGpu* nbnxmGpu)
315 {
316     if (const char* env = getenv("GMX_NB_MIN_CI"))
317     {
318         char* end = nullptr;
319
320         int minimumIlistCount = strtol(env, &end, 10);
321         if (!end || (*end != 0) || minimumIlistCount < 0)
322         {
323             gmx_fatal(
324                     FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
325         }
326
327         if (debug)
328         {
329             fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n", minimumIlistCount);
330         }
331         return minimumIlistCount;
332     }
333     else
334     {
335         int minimumIlistCount = gpu_min_ci_balanced(nbnxmGpu);
336         if (debug)
337         {
338             fprintf(debug,
339                     "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU "
340                     "multi-processors)\n",
341                     minimumIlistCount);
342         }
343         return minimumIlistCount;
344     }
345 }
346
347 static int getENbnxnInitCombRule(const t_forcerec& forcerec)
348 {
349     if (forcerec.ic->vdwtype == VanDerWaalsType::Cut
350         && (forcerec.ic->vdw_modifier == InteractionModifiers::None
351             || forcerec.ic->vdw_modifier == InteractionModifiers::PotShift)
352         && getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
353     {
354         /* Plain LJ cut-off: we can optimize with combination rules */
355         return enbnxninitcombruleDETECT;
356     }
357     else if (forcerec.ic->vdwtype == VanDerWaalsType::Pme)
358     {
359         /* LJ-PME: we need to use a combination rule for the grid */
360         if (forcerec.ljpme_combination_rule == LongRangeVdW::Geom)
361         {
362             return enbnxninitcombruleGEOM;
363         }
364         else
365         {
366             return enbnxninitcombruleLB;
367         }
368     }
369     else
370     {
371         /* We use a full combination matrix: no rule required */
372         return enbnxninitcombruleNONE;
373     }
374 }
375
376 std::unique_ptr<nonbonded_verlet_t> init_nb_verlet(const gmx::MDLogger& mdlog,
377                                                    const t_inputrec&    inputrec,
378                                                    const t_forcerec&    forcerec,
379                                                    const t_commrec*     commrec,
380                                                    const gmx_hw_info_t& hardwareInfo,
381                                                    bool                 useGpuForNonbonded,
382                                                    const gmx::DeviceStreamManager* deviceStreamManager,
383                                                    const gmx_mtop_t&               mtop,
384                                                    matrix                          box,
385                                                    gmx_wallcycle*                  wcycle)
386 {
387     const bool emulateGpu = (getenv("GMX_EMULATE_GPU") != nullptr);
388
389     GMX_RELEASE_ASSERT(!(emulateGpu && useGpuForNonbonded),
390                        "When GPU emulation is active, there cannot be a GPU assignment");
391
392     NonbondedResource nonbondedResource;
393     if (useGpuForNonbonded)
394     {
395         nonbondedResource = NonbondedResource::Gpu;
396     }
397     else if (emulateGpu)
398     {
399         nonbondedResource = NonbondedResource::EmulateGpu;
400     }
401     else
402     {
403         nonbondedResource = NonbondedResource::Cpu;
404     }
405
406     Nbnxm::KernelSetup kernelSetup = pick_nbnxn_kernel(
407             mdlog, forcerec.use_simd_kernels, hardwareInfo, nonbondedResource, inputrec);
408
409     const bool haveMultipleDomains = havePPDomainDecomposition(commrec);
410
411     bool bFEP_NonBonded = (forcerec.efep != FreeEnergyPerturbationType::No)
412                           && haveFepPerturbedNBInteractions(mtop);
413     PairlistParams pairlistParams(
414             kernelSetup.kernelType, bFEP_NonBonded, inputrec.rlist, haveMultipleDomains);
415
416     setupDynamicPairlistPruning(mdlog, inputrec, mtop, box, *forcerec.ic, &pairlistParams);
417
418     const int enbnxninitcombrule = getENbnxnInitCombRule(forcerec);
419
420     auto pinPolicy = (useGpuForNonbonded ? gmx::PinningPolicy::PinnedIfSupported
421                                          : gmx::PinningPolicy::CannotBePinned);
422
423     int mimimumNumEnergyGroupNonbonded = inputrec.opts.ngener;
424     if (inputrec.opts.ngener - inputrec.nwall == 1)
425     {
426         /* We have only one non-wall energy group, we do not need energy group
427          * support in the non-bondeds kernels, since all non-bonded energy
428          * contributions go to the first element of the energy group matrix.
429          */
430         mimimumNumEnergyGroupNonbonded = 1;
431     }
432
433     auto nbat = std::make_unique<nbnxn_atomdata_t>(
434             pinPolicy,
435             mdlog,
436             kernelSetup.kernelType,
437             enbnxninitcombrule,
438             forcerec.ntype,
439             forcerec.nbfp,
440             mimimumNumEnergyGroupNonbonded,
441             (useGpuForNonbonded || emulateGpu) ? 1 : gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded));
442
443     NbnxmGpu* gpu_nbv                          = nullptr;
444     int       minimumIlistCountForGpuBalancing = 0;
445     if (useGpuForNonbonded)
446     {
447         /* init the NxN GPU data; the last argument tells whether we'll have
448          * both local and non-local NB calculation on GPU */
449         GMX_RELEASE_ASSERT(
450                 (deviceStreamManager != nullptr),
451                 "Device stream manager should be initialized in order to use GPU for non-bonded.");
452         gpu_nbv = gpu_init(
453                 *deviceStreamManager, forcerec.ic.get(), pairlistParams, nbat.get(), haveMultipleDomains);
454
455         minimumIlistCountForGpuBalancing = getMinimumIlistCountForGpuBalancing(gpu_nbv);
456     }
457
458     auto pairlistSets = std::make_unique<PairlistSets>(
459             pairlistParams, haveMultipleDomains, minimumIlistCountForGpuBalancing);
460
461     auto pairSearch = std::make_unique<PairSearch>(
462             inputrec.pbcType,
463             EI_TPI(inputrec.eI),
464             haveDDAtomOrdering(*commrec) ? &commrec->dd->numCells : nullptr,
465             haveDDAtomOrdering(*commrec) ? domdec_zones(commrec->dd) : nullptr,
466             pairlistParams.pairlistType,
467             bFEP_NonBonded,
468             gmx_omp_nthreads_get(ModuleMultiThread::Pairsearch),
469             pinPolicy);
470
471     return std::make_unique<nonbonded_verlet_t>(
472             std::move(pairlistSets), std::move(pairSearch), std::move(nbat), kernelSetup, gpu_nbv, wcycle);
473 }
474
475 } // namespace Nbnxm
476
477 nonbonded_verlet_t::nonbonded_verlet_t(std::unique_ptr<PairlistSets>     pairlistSets,
478                                        std::unique_ptr<PairSearch>       pairSearch,
479                                        std::unique_ptr<nbnxn_atomdata_t> nbat_in,
480                                        const Nbnxm::KernelSetup&         kernelSetup,
481                                        NbnxmGpu*                         gpu_nbv_ptr,
482                                        gmx_wallcycle*                    wcycle) :
483     pairlistSets_(std::move(pairlistSets)),
484     pairSearch_(std::move(pairSearch)),
485     nbat(std::move(nbat_in)),
486     kernelSetup_(kernelSetup),
487     wcycle_(wcycle),
488     foreignEnergyGroups_(std::make_unique<gmx_grppairener_t>(nbat->params().nenergrp)),
489     gpu_nbv(gpu_nbv_ptr)
490 {
491     GMX_RELEASE_ASSERT(pairlistSets_, "Need valid pairlistSets");
492     GMX_RELEASE_ASSERT(pairSearch_, "Need valid search object");
493     GMX_RELEASE_ASSERT(nbat, "Need valid atomdata object");
494 }
495
496 nonbonded_verlet_t::~nonbonded_verlet_t()
497 {
498     Nbnxm::gpu_free(gpu_nbv);
499 }