Remove logging from hardware detection
[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/hardware/cpuinfo.h"
52 #include "gromacs/hardware/device_management.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/inmemoryserializer.h"
64 #include "gromacs/utility/logger.h"
65 #include "gromacs/utility/mutex.h"
66 #include "gromacs/utility/physicalnodecommunicator.h"
67
68 #include "architecture.h"
69 #include "device_information.h"
70
71 #ifdef HAVE_UNISTD_H
72 #    include <unistd.h> // sysconf()
73 #endif
74
75 gmx_hw_info_t::gmx_hw_info_t(std::unique_ptr<gmx::CpuInfo>          cpuInfo,
76                              std::unique_ptr<gmx::HardwareTopology> hardwareTopology) :
77     cpuInfo(std::move(cpuInfo)),
78     hardwareTopology(std::move(hardwareTopology))
79 {
80 }
81
82 gmx_hw_info_t::~gmx_hw_info_t() = default;
83
84 namespace gmx
85 {
86
87 //! Convenience macro to help us avoid ifdefs each time we use sysconf
88 #if !defined(_SC_NPROCESSORS_ONLN) && defined(_SC_NPROC_ONLN)
89 #    define _SC_NPROCESSORS_ONLN _SC_NPROC_ONLN
90 #endif
91
92 //! Convenience macro to help us avoid ifdefs each time we use sysconf
93 #if !defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROC_CONF)
94 #    define _SC_NPROCESSORS_CONF _SC_NPROC_CONF
95 #endif
96
97 /*! \brief The result of device detection
98  *
99  * Note that non-functional device detection still produces
100  * a detection result, ie. of no devices. This might not be
101  * what the user wanted, so it makes sense to log later when
102  * that is possible. */
103 struct DeviceDetectionResult
104 {
105     //! The device information detected
106     std::vector<std::unique_ptr<DeviceInformation>> deviceInfoList_;
107     //! Container of possible warnings to issue when that is possible
108     std::vector<std::string> deviceDetectionWarnings_;
109 };
110
111 /*! \brief Detect GPUs when that makes sense to attempt.
112  *
113  * \param[in]  physicalNodeComm  The communicator across this physical node
114  * \return The result of the detection, perhaps including diagnostic messages
115  *         to issue later.
116  *
117  * \todo Coordinating the efficient detection of devices across
118  * multiple ranks per node should be separated from the lower-level
119  * hardware detection. See
120  * https://gitlab.com/gromacs/gromacs/-/issues/3650.
121  */
122 static DeviceDetectionResult detectAllDeviceInformation(const PhysicalNodeCommunicator& physicalNodeComm)
123 {
124     DeviceDetectionResult deviceDetectionResult;
125
126     if (!isDeviceDetectionEnabled())
127     {
128         return deviceDetectionResult;
129     }
130
131     std::string errorMessage;
132
133     bool isMasterRankOfPhysicalNode = true;
134 #if GMX_LIB_MPI
135     isMasterRankOfPhysicalNode = (physicalNodeComm.rank_ == 0);
136 #else
137     // Without an MPI library, this process is trivially the only one
138     // on the physical node. This code runs before e.g. thread-MPI
139     // ranks are spawned, so detection is race-free by construction.
140     // Read-only access is enforced with providing those ranks with a
141     // handle to a const object, so usage is also free of races.
142     GMX_UNUSED_VALUE(physicalNodeComm);
143     isMasterRankOfPhysicalNode = true;
144 #endif
145
146     /* The SYCL and OpenCL support requires us to run detection on all
147      * ranks.
148      *
149      * With CUDA we don't need to, and prefer to detect on one rank
150      * and send the information to the other ranks over MPI. This
151      * avoids creating a start-up bottleneck with each MPI rank on a
152      * node making the same GPU API calls. */
153     constexpr bool allRanksMustDetectGpus = (GMX_GPU_OPENCL != 0 || GMX_GPU_SYCL != 0);
154     bool           gpusCanBeDetected      = false;
155     if (isMasterRankOfPhysicalNode || allRanksMustDetectGpus)
156     {
157         std::string errorMessage;
158         gpusCanBeDetected = isDeviceDetectionFunctional(&errorMessage);
159         if (!gpusCanBeDetected)
160         {
161             deviceDetectionResult.deviceDetectionWarnings_.emplace_back(
162                     "Detection of GPUs failed. The API reported:\n" + errorMessage);
163         }
164     }
165
166     if (gpusCanBeDetected)
167     {
168         deviceDetectionResult.deviceInfoList_ = findDevices();
169         // No need to tell the user anything at this point, they get a
170         // hardware report later.
171     }
172
173 #if GMX_LIB_MPI
174     if (!allRanksMustDetectGpus && (physicalNodeComm.size_ > 1))
175     {
176         // Master rank must serialize the device information list and
177         // send it to the other ranks on this node.
178         std::vector<char> buffer;
179         int               sizeOfBuffer;
180         if (isMasterRankOfPhysicalNode)
181         {
182             gmx::InMemorySerializer writer;
183             serializeDeviceInformations(deviceDetectionResult.deviceInfoList_, &writer);
184             buffer       = writer.finishAndGetBuffer();
185             sizeOfBuffer = buffer.size();
186         }
187         // Ensure all ranks agree on the size of list to be sent
188         MPI_Bcast(&sizeOfBuffer, 1, MPI_INT, 0, physicalNodeComm.comm_);
189         buffer.resize(sizeOfBuffer);
190         if (!buffer.empty())
191         {
192             // Send the list and deserialize it
193             MPI_Bcast(buffer.data(), buffer.size(), MPI_BYTE, 0, physicalNodeComm.comm_);
194             if (!isMasterRankOfPhysicalNode)
195             {
196                 gmx::InMemoryDeserializer reader(buffer, false);
197                 deviceDetectionResult.deviceInfoList_ = deserializeDeviceInformations(&reader);
198             }
199         }
200     }
201 #endif
202     return deviceDetectionResult;
203 }
204
205 //! Reduce the locally collected \p hardwareInfo over MPI ranks
206 static void gmx_collect_hardware_mpi(const gmx::CpuInfo&              cpuInfo,
207                                      const PhysicalNodeCommunicator&  physicalNodeComm,
208                                      compat::not_null<gmx_hw_info_t*> hardwareInfo)
209 {
210     const int ncore = hardwareInfo->hardwareTopology->numberOfCores();
211     /* Zen1 is assumed for:
212      * - family=23 with the below listed models;
213      * - Hygon as vendor.
214      */
215     const bool cpuIsAmdZen1 = ((cpuInfo.vendor() == CpuInfo::Vendor::Amd && cpuInfo.family() == 23
216                                 && (cpuInfo.model() == 1 || cpuInfo.model() == 17
217                                     || cpuInfo.model() == 8 || cpuInfo.model() == 24))
218                                || cpuInfo.vendor() == CpuInfo::Vendor::Hygon);
219
220     int numCompatibleDevices = getCompatibleDevices(hardwareInfo->deviceInfoList).size();
221 #if GMX_LIB_MPI
222     int nhwthread;
223     int gpu_hash;
224
225     nhwthread = hardwareInfo->nthreads_hw_avail;
226     /* Create a unique hash of the GPU type(s) in this node */
227     gpu_hash = 0;
228     /* Here it might be better to only loop over the compatible GPU, but we
229      * don't have that information available and it would also require
230      * removing the device ID from the device info string.
231      */
232     for (const auto& deviceInfo : hardwareInfo->deviceInfoList)
233     {
234         /* Since the device ID is incorporated in the hash, the order of
235          * the GPUs affects the hash. Also two identical GPUs won't give
236          * a gpu_hash of zero after XORing.
237          */
238         std::string deviceInfoString = getDeviceInformationString(*deviceInfo);
239         gpu_hash ^= gmx_string_fullhash_func(deviceInfoString.c_str(), gmx_string_hash_init);
240     }
241
242     constexpr int                      numElementsCounts = 4;
243     std::array<int, numElementsCounts> countsReduced;
244     {
245         std::array<int, numElementsCounts> countsLocal = { { 0 } };
246         // Organize to sum values from only one rank within each node,
247         // so we get the sum over all nodes.
248         bool isMasterRankOfPhysicalNode = (physicalNodeComm.rank_ == 0);
249         if (isMasterRankOfPhysicalNode)
250         {
251             countsLocal[0] = 1;
252             countsLocal[1] = ncore;
253             countsLocal[2] = nhwthread;
254             countsLocal[3] = numCompatibleDevices;
255         }
256
257         MPI_Allreduce(countsLocal.data(), countsReduced.data(), countsLocal.size(), MPI_INT,
258                       MPI_SUM, MPI_COMM_WORLD);
259     }
260
261     constexpr int                   numElementsMax = 11;
262     std::array<int, numElementsMax> maxMinReduced;
263     {
264         std::array<int, numElementsMax> maxMinLocal;
265         /* Store + and - values for all ranks,
266          * so we can get max+min with one MPI call.
267          */
268         maxMinLocal[0]  = ncore;
269         maxMinLocal[1]  = nhwthread;
270         maxMinLocal[2]  = numCompatibleDevices;
271         maxMinLocal[3]  = static_cast<int>(gmx::simdSuggested(cpuInfo));
272         maxMinLocal[4]  = gpu_hash;
273         maxMinLocal[5]  = -maxMinLocal[0];
274         maxMinLocal[6]  = -maxMinLocal[1];
275         maxMinLocal[7]  = -maxMinLocal[2];
276         maxMinLocal[8]  = -maxMinLocal[3];
277         maxMinLocal[9]  = -maxMinLocal[4];
278         maxMinLocal[10] = (cpuIsAmdZen1 ? 1 : 0);
279
280         MPI_Allreduce(maxMinLocal.data(), maxMinReduced.data(), maxMinLocal.size(), MPI_INT,
281                       MPI_MAX, MPI_COMM_WORLD);
282     }
283
284     hardwareInfo->nphysicalnode       = countsReduced[0];
285     hardwareInfo->ncore_tot           = countsReduced[1];
286     hardwareInfo->ncore_min           = -maxMinReduced[5];
287     hardwareInfo->ncore_max           = maxMinReduced[0];
288     hardwareInfo->nhwthread_tot       = countsReduced[2];
289     hardwareInfo->nhwthread_min       = -maxMinReduced[6];
290     hardwareInfo->nhwthread_max       = maxMinReduced[1];
291     hardwareInfo->ngpu_compatible_tot = countsReduced[3];
292     hardwareInfo->ngpu_compatible_min = -maxMinReduced[7];
293     hardwareInfo->ngpu_compatible_max = maxMinReduced[2];
294     hardwareInfo->simd_suggest_min    = -maxMinReduced[8];
295     hardwareInfo->simd_suggest_max    = maxMinReduced[3];
296     hardwareInfo->bIdenticalGPUs      = (maxMinReduced[4] == -maxMinReduced[9]);
297     hardwareInfo->haveAmdZen1Cpu      = (maxMinReduced[10] > 0);
298 #else
299     /* All ranks use the same pointer, protected by a mutex in the caller */
300     hardwareInfo->nphysicalnode       = 1;
301     hardwareInfo->ncore_tot           = ncore;
302     hardwareInfo->ncore_min           = ncore;
303     hardwareInfo->ncore_max           = ncore;
304     hardwareInfo->nhwthread_tot       = hardwareInfo->nthreads_hw_avail;
305     hardwareInfo->nhwthread_min       = hardwareInfo->nthreads_hw_avail;
306     hardwareInfo->nhwthread_max       = hardwareInfo->nthreads_hw_avail;
307     hardwareInfo->ngpu_compatible_tot = numCompatibleDevices;
308     hardwareInfo->ngpu_compatible_min = numCompatibleDevices;
309     hardwareInfo->ngpu_compatible_max = numCompatibleDevices;
310     hardwareInfo->simd_suggest_min    = static_cast<int>(simdSuggested(cpuInfo));
311     hardwareInfo->simd_suggest_max    = static_cast<int>(simdSuggested(cpuInfo));
312     hardwareInfo->bIdenticalGPUs      = TRUE;
313     hardwareInfo->haveAmdZen1Cpu      = cpuIsAmdZen1;
314     GMX_UNUSED_VALUE(physicalNodeComm);
315 #endif
316 }
317
318 /*! \brief Utility that does dummy computing for max 2 seconds to spin up cores
319  *
320  *  This routine will check the number of cores configured and online
321  *  (using sysconf), and the spins doing dummy compute operations for up to
322  *  2 seconds, or until all cores have come online. This can be used prior to
323  *  hardware detection for platforms that take unused processors offline.
324  *
325  *  This routine will not throw exceptions. In principle it should be
326  *  declared noexcept, but at least icc 19.1 and 21-beta08 with the
327  *  libstdc++-7.5 has difficulty implementing a std::vector of
328  *  std::thread started with this function when declared noexcept. It
329  *  is not clear whether the problem is the compiler or the standard
330  *  library. Fortunately, this function is not performance sensitive,
331  *  and only runs on platforms other than x86 and POWER (ie ARM),
332  *  so the possible overhead introduced by omitting noexcept is not
333  *  important.
334  */
335 static void spinUpCore()
336 {
337 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROCESSORS_ONLN)
338     float dummy           = 0.1;
339     int   countConfigured = sysconf(_SC_NPROCESSORS_CONF);    // noexcept
340     auto  start           = std::chrono::steady_clock::now(); // noexcept
341
342     while (sysconf(_SC_NPROCESSORS_ONLN) < countConfigured
343            && std::chrono::steady_clock::now() - start < std::chrono::seconds(2))
344     {
345         for (int i = 1; i < 10000; i++)
346         {
347             dummy /= i;
348         }
349     }
350
351     if (dummy < 0)
352     {
353         printf("This cannot happen, but prevents loop from being optimized away.");
354     }
355 #endif
356 }
357
358 /*! \brief Prepare the system before hardware topology detection
359  *
360  * This routine should perform any actions we want to put the system in a state
361  * where we want it to be before detecting the hardware topology. For most
362  * processors there is nothing to do, but some architectures (in particular ARM)
363  * have support for taking configured cores offline, which will make them disappear
364  * from the online processor count.
365  *
366  * This routine checks if there is a mismatch between the number of cores
367  * configured and online, and in that case we issue a small workload that
368  * attempts to wake sleeping cores before doing the actual detection.
369  *
370  * This type of mismatch can also occur for x86 or PowerPC on Linux, if SMT has only
371  * been disabled in the kernel (rather than bios). Since those cores will never
372  * come online automatically, we currently skip this test for x86 & PowerPC to
373  * avoid wasting 2 seconds. We also skip the test if there is no thread support.
374  *
375  * \note Cores will sleep relatively quickly again, so it's important to issue
376  *       the real detection code directly after this routine.
377  */
378 static void hardwareTopologyPrepareDetection()
379 {
380 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) \
381         && (defined(THREAD_PTHREADS) || defined(THREAD_WINDOWS))
382
383     // Modify this conditional when/if x86 or PowerPC starts to sleep some cores
384     if (c_architecture != Architecture::X86 && c_architecture != Architecture::PowerPC)
385     {
386         int                      countConfigured = sysconf(_SC_NPROCESSORS_CONF);
387         std::vector<std::thread> workThreads(countConfigured);
388
389         for (auto& t : workThreads)
390         {
391             t = std::thread(spinUpCore);
392         }
393
394         for (auto& t : workThreads)
395         {
396             t.join();
397         }
398     }
399 #endif
400 }
401
402 void hardwareTopologyDoubleCheckDetection(const gmx::MDLogger gmx_unused& mdlog,
403                                           const gmx::HardwareTopology gmx_unused& hardwareTopology)
404 {
405 #if defined HAVE_SYSCONF && defined(_SC_NPROCESSORS_CONF)
406     if (hardwareTopology.supportLevel() < gmx::HardwareTopology::SupportLevel::LogicalProcessorCount)
407     {
408         return;
409     }
410
411     int countFromDetection = hardwareTopology.machine().logicalProcessorCount;
412     int countConfigured    = sysconf(_SC_NPROCESSORS_CONF);
413
414     /* BIOS, kernel or user actions can take physical processors
415      * offline. We already cater for the some of the cases inside the hardwareToplogy
416      * by trying to spin up cores just before we detect, but there could be other
417      * cases where it is worthwhile to hint that there might be more resources available.
418      */
419     if (countConfigured >= 0 && countConfigured != countFromDetection)
420     {
421         GMX_LOG(mdlog.info)
422                 .appendTextFormatted(
423                         "Note: %d CPUs configured, but only %d were detected to be online.\n",
424                         countConfigured, countFromDetection);
425
426         if (c_architecture == Architecture::X86 && countConfigured == 2 * countFromDetection)
427         {
428             GMX_LOG(mdlog.info)
429                     .appendText(
430                             "      X86 Hyperthreading is likely disabled; enable it for better "
431                             "performance.");
432         }
433         // For PowerPC (likely Power8) it is possible to set SMT to either 2,4, or 8-way hardware threads.
434         // We only warn if it is completely disabled since default performance drops with SMT8.
435         if (c_architecture == Architecture::PowerPC && countConfigured == 8 * countFromDetection)
436         {
437             GMX_LOG(mdlog.info)
438                     .appendText(
439                             "      PowerPC SMT is likely disabled; enable SMT2/SMT4 for better "
440                             "performance.");
441         }
442     }
443 #else
444     GMX_UNUSED_VALUE(mdlog);
445     GMX_UNUSED_VALUE(hardwareTopology);
446 #endif
447 }
448
449 std::unique_ptr<gmx_hw_info_t> gmx_detect_hardware(const PhysicalNodeCommunicator& physicalNodeComm)
450 {
451     // Make the new hardwareInfo in a temporary.
452     hardwareTopologyPrepareDetection();
453
454     // TODO: We should also do CPU hardware detection only once on each
455     // physical node and broadcast it, instead of doing it on every MPI rank.
456     auto hardwareInfo = std::make_unique<gmx_hw_info_t>(
457             std::make_unique<CpuInfo>(CpuInfo::detect()),
458             std::make_unique<HardwareTopology>(HardwareTopology::detect()));
459
460     // TODO: Get rid of this altogether.
461     hardwareInfo->nthreads_hw_avail = hardwareInfo->hardwareTopology->machine().logicalProcessorCount;
462
463     // Detect GPUs
464     // Open a nested scope so no temporary variables can
465     // be mis-used later.
466     {
467         DeviceDetectionResult deviceDetectionResult = detectAllDeviceInformation(physicalNodeComm);
468         hardwareInfo->deviceInfoList.swap(deviceDetectionResult.deviceInfoList_);
469         std::swap(hardwareInfo->hardwareDetectionWarnings_, deviceDetectionResult.deviceDetectionWarnings_);
470     }
471
472     gmx_collect_hardware_mpi(*hardwareInfo->cpuInfo, physicalNodeComm, compat::make_not_null(hardwareInfo));
473
474     return hardwareInfo;
475 }
476
477 void logHardwareDetectionWarnings(const gmx::MDLogger& mdlog, const gmx_hw_info_t& hardwareInformation)
478 {
479     for (const std::string& warningString : hardwareInformation.hardwareDetectionWarnings_)
480     {
481         GMX_LOG(mdlog.warning).asParagraph().appendText(warningString);
482     }
483 }
484
485 } // namespace gmx