e88d642dc328f334139236e69dbe0f7c01d0814c
[alexxy/gromacs.git] / src / gromacs / hardware / detecthardware.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 #include "gmxpre.h"
37
38 #include "detecthardware.h"
39
40 #include "config.h"
41
42 #include <algorithm>
43 #include <array>
44 #include <chrono>
45 #include <memory>
46 #include <string>
47 #include <thread>
48 #include <vector>
49
50 #include "gromacs/compat/pointers.h"
51 #include "gromacs/gpu_utils/gpu_utils.h"
52 #include "gromacs/hardware/cpuinfo.h"
53 #include "gromacs/hardware/hardwaretopology.h"
54 #include "gromacs/hardware/hw_info.h"
55 #include "gromacs/simd/support.h"
56 #include "gromacs/utility/basedefinitions.h"
57 #include "gromacs/utility/basenetwork.h"
58 #include "gromacs/utility/baseversion.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/exceptions.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/gmxassert.h"
63 #include "gromacs/utility/logger.h"
64 #include "gromacs/utility/mutex.h"
65 #include "gromacs/utility/physicalnodecommunicator.h"
66
67 #include "architecture.h"
68
69 #ifdef HAVE_UNISTD_H
70 #    include <unistd.h> // sysconf()
71 #endif
72
73 gmx_hw_info_t::gmx_hw_info_t(std::unique_ptr<gmx::CpuInfo>          cpuInfo,
74                              std::unique_ptr<gmx::HardwareTopology> hardwareTopology) :
75     cpuInfo(std::move(cpuInfo)),
76     hardwareTopology(std::move(hardwareTopology))
77 {
78 }
79
80 gmx_hw_info_t::~gmx_hw_info_t()
81 {
82     free_gpu_info(&gpu_info);
83 }
84
85 namespace gmx
86 {
87
88 //! Convenience macro to help us avoid ifdefs each time we use sysconf
89 #if !defined(_SC_NPROCESSORS_ONLN) && defined(_SC_NPROC_ONLN)
90 #    define _SC_NPROCESSORS_ONLN _SC_NPROC_ONLN
91 #endif
92
93 //! Convenience macro to help us avoid ifdefs each time we use sysconf
94 #if !defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROC_CONF)
95 #    define _SC_NPROCESSORS_CONF _SC_NPROC_CONF
96 #endif
97
98 /*! \brief Information about the hardware of all nodes (common to all threads in this process).
99  *
100  * This information is constructed only when required, but thereafter
101  * its lifetime is that of the whole process, potentially across
102  * multiple successive simulation parts. It's wise to ensure that only
103  * one thread can create the information, but thereafter they can all
104  * read it without e.g. needing a std::shared_ptr to ensure its
105  * lifetime exceeds that of the thread. */
106 static std::unique_ptr<gmx_hw_info_t> g_hardwareInfo;
107 //! A mutex to protect the hwinfo structure
108 static Mutex g_hardwareInfoMutex;
109
110 //! Detect GPUs, if that makes sense to attempt.
111 static void gmx_detect_gpus(const gmx::MDLogger&             mdlog,
112                             const PhysicalNodeCommunicator&  physicalNodeComm,
113                             compat::not_null<gmx_hw_info_t*> hardwareInfo)
114 {
115     hardwareInfo->gpu_info.bDetectGPUs = canPerformGpuDetection();
116
117     if (!hardwareInfo->gpu_info.bDetectGPUs)
118     {
119         return;
120     }
121
122     bool isMasterRankOfPhysicalNode = true;
123 #if GMX_LIB_MPI
124     isMasterRankOfPhysicalNode = (physicalNodeComm.rank_ == 0);
125 #else
126     // We choose to run the detection only once with thread-MPI and
127     // use a mutex to enforce it.
128     GMX_UNUSED_VALUE(physicalNodeComm);
129     isMasterRankOfPhysicalNode = true;
130 #endif
131
132     /* The OpenCL support requires us to run detection on all ranks.
133      * With CUDA we don't need to, and prefer to detect on one rank
134      * and send the information to the other ranks over MPI. */
135     constexpr bool allRanksMustDetectGpus = (GMX_GPU_OPENCL != 0);
136     bool           gpusCanBeDetected      = false;
137     if (isMasterRankOfPhysicalNode || allRanksMustDetectGpus)
138     {
139         std::string errorMessage;
140         gpusCanBeDetected = isGpuDetectionFunctional(&errorMessage);
141         if (!gpusCanBeDetected)
142         {
143             GMX_LOG(mdlog.info)
144                     .asParagraph()
145                     .appendTextFormatted(
146                             "NOTE: Detection of GPUs failed. The API reported:\n"
147                             "      %s\n"
148                             "      GROMACS cannot run tasks on a GPU.",
149                             errorMessage.c_str());
150         }
151     }
152
153     if (gpusCanBeDetected)
154     {
155         findGpus(&hardwareInfo->gpu_info);
156         // No need to tell the user anything at this point, they get a
157         // hardware report later.
158     }
159
160 #if GMX_LIB_MPI
161     if (!allRanksMustDetectGpus)
162     {
163         /* Broadcast the GPU info to the other ranks within this node */
164         MPI_Bcast(&hardwareInfo->gpu_info.n_dev, 1, MPI_INT, 0, physicalNodeComm.comm_);
165
166         if (hardwareInfo->gpu_info.n_dev > 0)
167         {
168             int dev_size;
169
170             dev_size = hardwareInfo->gpu_info.n_dev * sizeof_gpu_dev_info();
171
172             if (!isMasterRankOfPhysicalNode)
173             {
174                 hardwareInfo->gpu_info.deviceInfo = (struct DeviceInformation*)malloc(dev_size);
175             }
176             MPI_Bcast(hardwareInfo->gpu_info.deviceInfo, dev_size, MPI_BYTE, 0, physicalNodeComm.comm_);
177             MPI_Bcast(&hardwareInfo->gpu_info.n_dev_compatible, 1, MPI_INT, 0, physicalNodeComm.comm_);
178         }
179     }
180 #endif
181 }
182
183 //! Reduce the locally collected \p hardwareInfo over MPI ranks
184 static void gmx_collect_hardware_mpi(const gmx::CpuInfo&              cpuInfo,
185                                      const PhysicalNodeCommunicator&  physicalNodeComm,
186                                      compat::not_null<gmx_hw_info_t*> hardwareInfo)
187 {
188     const int ncore = hardwareInfo->hardwareTopology->numberOfCores();
189     /* Zen1 is assumed for:
190      * - family=23 with the below listed models;
191      * - Hygon as vendor.
192      */
193     const bool cpuIsAmdZen1 = ((cpuInfo.vendor() == CpuInfo::Vendor::Amd && cpuInfo.family() == 23
194                                 && (cpuInfo.model() == 1 || cpuInfo.model() == 17
195                                     || cpuInfo.model() == 8 || cpuInfo.model() == 24))
196                                || cpuInfo.vendor() == CpuInfo::Vendor::Hygon);
197 #if GMX_LIB_MPI
198     int nhwthread, ngpu, i;
199     int gpu_hash;
200
201     nhwthread = hardwareInfo->nthreads_hw_avail;
202     ngpu      = hardwareInfo->gpu_info.n_dev_compatible;
203     /* Create a unique hash of the GPU type(s) in this node */
204     gpu_hash = 0;
205     /* Here it might be better to only loop over the compatible GPU, but we
206      * don't have that information available and it would also require
207      * removing the device ID from the device info string.
208      */
209     for (i = 0; i < hardwareInfo->gpu_info.n_dev; i++)
210     {
211         char stmp[STRLEN];
212
213         /* Since the device ID is incorporated in the hash, the order of
214          * the GPUs affects the hash. Also two identical GPUs won't give
215          * a gpu_hash of zero after XORing.
216          */
217         get_gpu_device_info_string(stmp, hardwareInfo->gpu_info, i);
218         gpu_hash ^= gmx_string_fullhash_func(stmp, gmx_string_hash_init);
219     }
220
221     constexpr int                      numElementsCounts = 4;
222     std::array<int, numElementsCounts> countsReduced;
223     {
224         std::array<int, numElementsCounts> countsLocal = { { 0 } };
225         // Organize to sum values from only one rank within each node,
226         // so we get the sum over all nodes.
227         bool isMasterRankOfPhysicalNode = (physicalNodeComm.rank_ == 0);
228         if (isMasterRankOfPhysicalNode)
229         {
230             countsLocal[0] = 1;
231             countsLocal[1] = ncore;
232             countsLocal[2] = nhwthread;
233             countsLocal[3] = ngpu;
234         }
235
236         MPI_Allreduce(countsLocal.data(), countsReduced.data(), countsLocal.size(), MPI_INT,
237                       MPI_SUM, MPI_COMM_WORLD);
238     }
239
240     constexpr int                   numElementsMax = 11;
241     std::array<int, numElementsMax> maxMinReduced;
242     {
243         std::array<int, numElementsMax> maxMinLocal;
244         /* Store + and - values for all ranks,
245          * so we can get max+min with one MPI call.
246          */
247         maxMinLocal[0]  = ncore;
248         maxMinLocal[1]  = nhwthread;
249         maxMinLocal[2]  = ngpu;
250         maxMinLocal[3]  = static_cast<int>(gmx::simdSuggested(cpuInfo));
251         maxMinLocal[4]  = gpu_hash;
252         maxMinLocal[5]  = -maxMinLocal[0];
253         maxMinLocal[6]  = -maxMinLocal[1];
254         maxMinLocal[7]  = -maxMinLocal[2];
255         maxMinLocal[8]  = -maxMinLocal[3];
256         maxMinLocal[9]  = -maxMinLocal[4];
257         maxMinLocal[10] = (cpuIsAmdZen1 ? 1 : 0);
258
259         MPI_Allreduce(maxMinLocal.data(), maxMinReduced.data(), maxMinLocal.size(), MPI_INT,
260                       MPI_MAX, MPI_COMM_WORLD);
261     }
262
263     hardwareInfo->nphysicalnode       = countsReduced[0];
264     hardwareInfo->ncore_tot           = countsReduced[1];
265     hardwareInfo->ncore_min           = -maxMinReduced[5];
266     hardwareInfo->ncore_max           = maxMinReduced[0];
267     hardwareInfo->nhwthread_tot       = countsReduced[2];
268     hardwareInfo->nhwthread_min       = -maxMinReduced[6];
269     hardwareInfo->nhwthread_max       = maxMinReduced[1];
270     hardwareInfo->ngpu_compatible_tot = countsReduced[3];
271     hardwareInfo->ngpu_compatible_min = -maxMinReduced[7];
272     hardwareInfo->ngpu_compatible_max = maxMinReduced[2];
273     hardwareInfo->simd_suggest_min    = -maxMinReduced[8];
274     hardwareInfo->simd_suggest_max    = maxMinReduced[3];
275     hardwareInfo->bIdenticalGPUs      = (maxMinReduced[4] == -maxMinReduced[9]);
276     hardwareInfo->haveAmdZen1Cpu      = (maxMinReduced[10] > 0);
277 #else
278     /* All ranks use the same pointer, protected by a mutex in the caller */
279     hardwareInfo->nphysicalnode       = 1;
280     hardwareInfo->ncore_tot           = ncore;
281     hardwareInfo->ncore_min           = ncore;
282     hardwareInfo->ncore_max           = ncore;
283     hardwareInfo->nhwthread_tot       = hardwareInfo->nthreads_hw_avail;
284     hardwareInfo->nhwthread_min       = hardwareInfo->nthreads_hw_avail;
285     hardwareInfo->nhwthread_max       = hardwareInfo->nthreads_hw_avail;
286     hardwareInfo->ngpu_compatible_tot = hardwareInfo->gpu_info.n_dev_compatible;
287     hardwareInfo->ngpu_compatible_min = hardwareInfo->gpu_info.n_dev_compatible;
288     hardwareInfo->ngpu_compatible_max = hardwareInfo->gpu_info.n_dev_compatible;
289     hardwareInfo->simd_suggest_min    = static_cast<int>(simdSuggested(cpuInfo));
290     hardwareInfo->simd_suggest_max    = static_cast<int>(simdSuggested(cpuInfo));
291     hardwareInfo->bIdenticalGPUs      = TRUE;
292     hardwareInfo->haveAmdZen1Cpu      = cpuIsAmdZen1;
293     GMX_UNUSED_VALUE(physicalNodeComm);
294 #endif
295 }
296
297 /*! \brief Utility that does dummy computing for max 2 seconds to spin up cores
298  *
299  *  This routine will check the number of cores configured and online
300  *  (using sysconf), and the spins doing dummy compute operations for up to
301  *  2 seconds, or until all cores have come online. This can be used prior to
302  *  hardware detection for platforms that take unused processors offline.
303  *
304  *  This routine will not throw exceptions.
305  */
306 static void spinUpCore() noexcept
307 {
308 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROCESSORS_ONLN)
309     float dummy           = 0.1;
310     int   countConfigured = sysconf(_SC_NPROCESSORS_CONF);    // noexcept
311     auto  start           = std::chrono::steady_clock::now(); // noexcept
312
313     while (sysconf(_SC_NPROCESSORS_ONLN) < countConfigured
314            && std::chrono::steady_clock::now() - start < std::chrono::seconds(2))
315     {
316         for (int i = 1; i < 10000; i++)
317         {
318             dummy /= i;
319         }
320     }
321
322     if (dummy < 0)
323     {
324         printf("This cannot happen, but prevents loop from being optimized away.");
325     }
326 #endif
327 }
328
329 /*! \brief Prepare the system before hardware topology detection
330  *
331  * This routine should perform any actions we want to put the system in a state
332  * where we want it to be before detecting the hardware topology. For most
333  * processors there is nothing to do, but some architectures (in particular ARM)
334  * have support for taking configured cores offline, which will make them disappear
335  * from the online processor count.
336  *
337  * This routine checks if there is a mismatch between the number of cores
338  * configured and online, and in that case we issue a small workload that
339  * attempts to wake sleeping cores before doing the actual detection.
340  *
341  * This type of mismatch can also occur for x86 or PowerPC on Linux, if SMT has only
342  * been disabled in the kernel (rather than bios). Since those cores will never
343  * come online automatically, we currently skip this test for x86 & PowerPC to
344  * avoid wasting 2 seconds. We also skip the test if there is no thread support.
345  *
346  * \note Cores will sleep relatively quickly again, so it's important to issue
347  *       the real detection code directly after this routine.
348  */
349 static void hardwareTopologyPrepareDetection()
350 {
351 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) \
352         && (defined(THREAD_PTHREADS) || defined(THREAD_WINDOWS))
353
354     // Modify this conditional when/if x86 or PowerPC starts to sleep some cores
355     if (c_architecture != Architecture::X86 && c_architecture != Architecture::PowerPC)
356     {
357         int                      countConfigured = sysconf(_SC_NPROCESSORS_CONF);
358         std::vector<std::thread> workThreads(countConfigured);
359
360         for (auto& t : workThreads)
361         {
362             t = std::thread(spinUpCore);
363         }
364
365         for (auto& t : workThreads)
366         {
367             t.join();
368         }
369     }
370 #endif
371 }
372
373 /*! \brief Sanity check hardware topology and print some notes to log
374  *
375  *  \param mdlog            Logger.
376  *  \param hardwareTopology Reference to hardwareTopology object.
377  */
378 static void hardwareTopologyDoubleCheckDetection(const gmx::MDLogger gmx_unused& mdlog,
379                                                  const gmx::HardwareTopology gmx_unused& hardwareTopology)
380 {
381 #if defined HAVE_SYSCONF && defined(_SC_NPROCESSORS_CONF)
382     if (hardwareTopology.supportLevel() < gmx::HardwareTopology::SupportLevel::LogicalProcessorCount)
383     {
384         return;
385     }
386
387     int countFromDetection = hardwareTopology.machine().logicalProcessorCount;
388     int countConfigured    = sysconf(_SC_NPROCESSORS_CONF);
389
390     /* BIOS, kernel or user actions can take physical processors
391      * offline. We already cater for the some of the cases inside the hardwareToplogy
392      * by trying to spin up cores just before we detect, but there could be other
393      * cases where it is worthwhile to hint that there might be more resources available.
394      */
395     if (countConfigured >= 0 && countConfigured != countFromDetection)
396     {
397         GMX_LOG(mdlog.info)
398                 .appendTextFormatted(
399                         "Note: %d CPUs configured, but only %d were detected to be online.\n",
400                         countConfigured, countFromDetection);
401
402         if (c_architecture == Architecture::X86 && countConfigured == 2 * countFromDetection)
403         {
404             GMX_LOG(mdlog.info)
405                     .appendText(
406                             "      X86 Hyperthreading is likely disabled; enable it for better "
407                             "performance.");
408         }
409         // For PowerPC (likely Power8) it is possible to set SMT to either 2,4, or 8-way hardware threads.
410         // We only warn if it is completely disabled since default performance drops with SMT8.
411         if (c_architecture == Architecture::PowerPC && countConfigured == 8 * countFromDetection)
412         {
413             GMX_LOG(mdlog.info)
414                     .appendText(
415                             "      PowerPC SMT is likely disabled; enable SMT2/SMT4 for better "
416                             "performance.");
417         }
418     }
419 #endif
420 }
421
422 gmx_hw_info_t* gmx_detect_hardware(const gmx::MDLogger& mdlog, const PhysicalNodeCommunicator& physicalNodeComm)
423 {
424     // By construction, only one thread ever runs hardware detection,
425     // but we may as well prevent issues arising if that would change.
426     // Taking the lock early ensures that exactly one thread can
427     // attempt to construct g_hardwareInfo.
428     lock_guard<Mutex> lock(g_hardwareInfoMutex);
429
430     // If we already have the information, just return a handle to it.
431     if (g_hardwareInfo != nullptr)
432     {
433         return g_hardwareInfo.get();
434     }
435
436     // Make the new hardwareInfo in a temporary.
437     hardwareTopologyPrepareDetection();
438
439     // TODO: We should also do CPU hardware detection only once on each
440     // physical node and broadcast it, instead of doing it on every MPI rank.
441     auto hardwareInfo = std::make_unique<gmx_hw_info_t>(
442             std::make_unique<CpuInfo>(CpuInfo::detect()),
443             std::make_unique<HardwareTopology>(HardwareTopology::detect()));
444
445     // If we detected the topology on this system, double-check that it makes sense
446     if (hardwareInfo->hardwareTopology->isThisSystem())
447     {
448         hardwareTopologyDoubleCheckDetection(mdlog, *hardwareInfo->hardwareTopology);
449     }
450
451     // TODO: Get rid of this altogether.
452     hardwareInfo->nthreads_hw_avail = hardwareInfo->hardwareTopology->machine().logicalProcessorCount;
453
454     // Detect GPUs
455     hardwareInfo->gpu_info.n_dev            = 0;
456     hardwareInfo->gpu_info.n_dev_compatible = 0;
457     hardwareInfo->gpu_info.deviceInfo       = nullptr;
458
459     gmx_detect_gpus(mdlog, physicalNodeComm, compat::make_not_null(hardwareInfo));
460     gmx_collect_hardware_mpi(*hardwareInfo->cpuInfo, physicalNodeComm, compat::make_not_null(hardwareInfo));
461
462     // Now that the temporary is fully constructed, swap it to become
463     // the real thing.
464     g_hardwareInfo.swap(hardwareInfo);
465
466     return g_hardwareInfo.get();
467 }
468
469 bool compatibleGpusFound(const gmx_gpu_info_t& gpu_info)
470 {
471     return gpu_info.n_dev_compatible > 0;
472 }
473
474 } // namespace gmx