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