Update mdrun message to reflect that Cuda CC>=3.5 is supported.
[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, The GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020,2021, 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 <memory>
45 #include <string>
46 #include <vector>
47
48 #include "gromacs/hardware/cpuinfo.h"
49 #include "gromacs/hardware/device_management.h"
50 #include "gromacs/hardware/hardwaretopology.h"
51 #include "gromacs/hardware/hw_info.h"
52 #include "gromacs/simd/support.h"
53 #include "gromacs/utility/basedefinitions.h"
54 #include "gromacs/utility/basenetwork.h"
55 #include "gromacs/utility/baseversion.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/inmemoryserializer.h"
61 #include "gromacs/utility/logger.h"
62 #include "gromacs/utility/physicalnodecommunicator.h"
63
64 #include "architecture.h"
65 #include "device_information.h"
66 #include "prepare_detection.h"
67
68 #ifdef HAVE_UNISTD_H
69 #    include <unistd.h> // sysconf()
70 #endif
71
72 gmx_hw_info_t::gmx_hw_info_t(std::unique_ptr<gmx::CpuInfo>          cpuInfo,
73                              std::unique_ptr<gmx::HardwareTopology> hardwareTopology) :
74     cpuInfo(std::move(cpuInfo)), hardwareTopology(std::move(hardwareTopology))
75 {
76 }
77
78 gmx_hw_info_t::~gmx_hw_info_t() = default;
79
80 namespace gmx
81 {
82
83 //! Convenience macro to help us avoid ifdefs each time we use sysconf
84 #if !defined(_SC_NPROCESSORS_ONLN) && defined(_SC_NPROC_ONLN)
85 #    define _SC_NPROCESSORS_ONLN _SC_NPROC_ONLN
86 #endif
87
88 //! Convenience macro to help us avoid ifdefs each time we use sysconf
89 #if !defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROC_CONF)
90 #    define _SC_NPROCESSORS_CONF _SC_NPROC_CONF
91 #endif
92
93 /*! \brief The result of device detection
94  *
95  * Note that non-functional device detection still produces
96  * a detection result, ie. of no devices. This might not be
97  * what the user wanted, so it makes sense to log later when
98  * that is possible. */
99 struct DeviceDetectionResult
100 {
101     //! The device information detected
102     std::vector<std::unique_ptr<DeviceInformation>> deviceInfoList_;
103     //! Container of possible warnings to issue when that is possible
104     std::vector<std::string> deviceDetectionWarnings_;
105 };
106
107 /*! \brief Detect GPUs when that makes sense to attempt.
108  *
109  * \param[in]  physicalNodeComm  The communicator across this physical node
110  * \return The result of the detection, perhaps including diagnostic messages
111  *         to issue later.
112  *
113  * \todo Coordinating the efficient detection of devices across
114  * multiple ranks per node should be separated from the lower-level
115  * hardware detection. See
116  * https://gitlab.com/gromacs/gromacs/-/issues/3650.
117  */
118 static DeviceDetectionResult detectAllDeviceInformation(const PhysicalNodeCommunicator& physicalNodeComm)
119 {
120     DeviceDetectionResult deviceDetectionResult;
121
122     if (!isDeviceDetectionEnabled())
123     {
124         return deviceDetectionResult;
125     }
126
127     std::string errorMessage;
128
129     bool isMasterRankOfPhysicalNode = true;
130 #if GMX_LIB_MPI
131     isMasterRankOfPhysicalNode = (physicalNodeComm.rank_ == 0);
132 #else
133     // Without an MPI library, this process is trivially the only one
134     // on the physical node. This code runs before e.g. thread-MPI
135     // ranks are spawned, so detection is race-free by construction.
136     // Read-only access is enforced with providing those ranks with a
137     // handle to a const object, so usage is also free of races.
138     GMX_UNUSED_VALUE(physicalNodeComm);
139     isMasterRankOfPhysicalNode        = true;
140 #endif
141
142     /* The SYCL and OpenCL support requires us to run detection on all
143      * ranks.
144      *
145      * With CUDA we don't need to, and prefer to detect on one rank
146      * and send the information to the other ranks over MPI. This
147      * avoids creating a start-up bottleneck with each MPI rank on a
148      * node making the same GPU API calls. */
149     constexpr bool allRanksMustDetectGpus = (GMX_GPU_OPENCL != 0 || GMX_GPU_SYCL != 0);
150     bool           gpusCanBeDetected      = false;
151     if (isMasterRankOfPhysicalNode || allRanksMustDetectGpus)
152     {
153         std::string errorMessage;
154         gpusCanBeDetected = isDeviceDetectionFunctional(&errorMessage);
155         if (!gpusCanBeDetected)
156         {
157             deviceDetectionResult.deviceDetectionWarnings_.emplace_back(
158                     "Detection of GPUs failed. The API reported:\n" + errorMessage);
159         }
160     }
161
162     if (gpusCanBeDetected)
163     {
164         deviceDetectionResult.deviceInfoList_ = findDevices();
165         // No need to tell the user anything at this point, they get a
166         // hardware report later.
167     }
168
169 #if GMX_LIB_MPI
170     if (!allRanksMustDetectGpus && (physicalNodeComm.size_ > 1))
171     {
172         // Master rank must serialize the device information list and
173         // send it to the other ranks on this node.
174         std::vector<char> buffer;
175         int               sizeOfBuffer;
176         if (isMasterRankOfPhysicalNode)
177         {
178             gmx::InMemorySerializer writer;
179             serializeDeviceInformations(deviceDetectionResult.deviceInfoList_, &writer);
180             buffer       = writer.finishAndGetBuffer();
181             sizeOfBuffer = buffer.size();
182         }
183         // Ensure all ranks agree on the size of list to be sent
184         MPI_Bcast(&sizeOfBuffer, 1, MPI_INT, 0, physicalNodeComm.comm_);
185         buffer.resize(sizeOfBuffer);
186         if (!buffer.empty())
187         {
188             // Send the list and deserialize it
189             MPI_Bcast(buffer.data(), buffer.size(), MPI_BYTE, 0, physicalNodeComm.comm_);
190             if (!isMasterRankOfPhysicalNode)
191             {
192                 gmx::InMemoryDeserializer reader(buffer, false);
193                 deviceDetectionResult.deviceInfoList_ = deserializeDeviceInformations(&reader);
194             }
195         }
196     }
197 #endif
198     return deviceDetectionResult;
199 }
200
201 //! Reduce the locally collected \p hardwareInfo over MPI ranks
202 static void gmx_collect_hardware_mpi(const gmx::CpuInfo&             cpuInfo,
203                                      const PhysicalNodeCommunicator& physicalNodeComm,
204                                      gmx_hw_info_t*                  hardwareInfo)
205 {
206     const int ncore = hardwareInfo->hardwareTopology->numberOfCores();
207     /* Zen1 is assumed for:
208      * - family=23 with the below listed models;
209      * - Hygon as vendor.
210      */
211     const bool cpuIsAmdZen1 = gmx::cpuIsAmdZen1(cpuInfo);
212
213     int numCompatibleDevices = getCompatibleDevices(hardwareInfo->deviceInfoList).size();
214 #if GMX_LIB_MPI
215     int nhwthread;
216     int gpu_hash;
217
218     nhwthread = hardwareInfo->nthreads_hw_avail;
219     /* Create a unique hash of the GPU type(s) in this node */
220     gpu_hash = 0;
221     /* Here it might be better to only loop over the compatible GPU, but we
222      * don't have that information available and it would also require
223      * removing the device ID from the device info string.
224      */
225     for (const auto& deviceInfo : hardwareInfo->deviceInfoList)
226     {
227         /* Since the device ID is incorporated in the hash, the order of
228          * the GPUs affects the hash. Also two identical GPUs won't give
229          * a gpu_hash of zero after XORing.
230          */
231         std::string deviceInfoString = getDeviceInformationString(*deviceInfo);
232         gpu_hash ^= gmx_string_fullhash_func(deviceInfoString.c_str(), gmx_string_hash_init);
233     }
234
235     constexpr int                      numElementsCounts = 4;
236     std::array<int, numElementsCounts> countsReduced;
237     {
238         std::array<int, numElementsCounts> countsLocal = { { 0 } };
239         // Organize to sum values from only one rank within each node,
240         // so we get the sum over all nodes.
241         bool isMasterRankOfPhysicalNode = (physicalNodeComm.rank_ == 0);
242         if (isMasterRankOfPhysicalNode)
243         {
244             countsLocal[0] = 1;
245             countsLocal[1] = ncore;
246             countsLocal[2] = nhwthread;
247             countsLocal[3] = numCompatibleDevices;
248         }
249
250         MPI_Allreduce(countsLocal.data(), countsReduced.data(), countsLocal.size(), MPI_INT, MPI_SUM, MPI_COMM_WORLD);
251     }
252
253     constexpr int                   numElementsMax = 11;
254     std::array<int, numElementsMax> maxMinReduced;
255     {
256         std::array<int, numElementsMax> maxMinLocal;
257         /* Store + and - values for all ranks,
258          * so we can get max+min with one MPI call.
259          */
260         maxMinLocal[0]  = ncore;
261         maxMinLocal[1]  = nhwthread;
262         maxMinLocal[2]  = numCompatibleDevices;
263         maxMinLocal[3]  = static_cast<int>(gmx::simdSuggested(cpuInfo));
264         maxMinLocal[4]  = gpu_hash;
265         maxMinLocal[5]  = -maxMinLocal[0];
266         maxMinLocal[6]  = -maxMinLocal[1];
267         maxMinLocal[7]  = -maxMinLocal[2];
268         maxMinLocal[8]  = -maxMinLocal[3];
269         maxMinLocal[9]  = -maxMinLocal[4];
270         maxMinLocal[10] = (cpuIsAmdZen1 ? 1 : 0);
271
272         MPI_Allreduce(maxMinLocal.data(), maxMinReduced.data(), maxMinLocal.size(), MPI_INT, MPI_MAX, MPI_COMM_WORLD);
273     }
274
275     hardwareInfo->nphysicalnode       = countsReduced[0];
276     hardwareInfo->ncore_tot           = countsReduced[1];
277     hardwareInfo->ncore_min           = -maxMinReduced[5];
278     hardwareInfo->ncore_max           = maxMinReduced[0];
279     hardwareInfo->nhwthread_tot       = countsReduced[2];
280     hardwareInfo->nhwthread_min       = -maxMinReduced[6];
281     hardwareInfo->nhwthread_max       = maxMinReduced[1];
282     hardwareInfo->ngpu_compatible_tot = countsReduced[3];
283     hardwareInfo->ngpu_compatible_min = -maxMinReduced[7];
284     hardwareInfo->ngpu_compatible_max = maxMinReduced[2];
285     hardwareInfo->simd_suggest_min    = -maxMinReduced[8];
286     hardwareInfo->simd_suggest_max    = maxMinReduced[3];
287     hardwareInfo->bIdenticalGPUs      = (maxMinReduced[4] == -maxMinReduced[9]);
288     hardwareInfo->haveAmdZen1Cpu      = (maxMinReduced[10] > 0);
289 #else
290     hardwareInfo->nphysicalnode       = 1;
291     hardwareInfo->ncore_tot           = ncore;
292     hardwareInfo->ncore_min           = ncore;
293     hardwareInfo->ncore_max           = ncore;
294     hardwareInfo->nhwthread_tot       = hardwareInfo->nthreads_hw_avail;
295     hardwareInfo->nhwthread_min       = hardwareInfo->nthreads_hw_avail;
296     hardwareInfo->nhwthread_max       = hardwareInfo->nthreads_hw_avail;
297     hardwareInfo->ngpu_compatible_tot = numCompatibleDevices;
298     hardwareInfo->ngpu_compatible_min = numCompatibleDevices;
299     hardwareInfo->ngpu_compatible_max = numCompatibleDevices;
300     hardwareInfo->simd_suggest_min    = static_cast<int>(simdSuggested(cpuInfo));
301     hardwareInfo->simd_suggest_max    = static_cast<int>(simdSuggested(cpuInfo));
302     hardwareInfo->bIdenticalGPUs      = TRUE;
303     hardwareInfo->haveAmdZen1Cpu      = cpuIsAmdZen1;
304     GMX_UNUSED_VALUE(physicalNodeComm);
305 #endif
306 }
307
308 void hardwareTopologyDoubleCheckDetection(const gmx::MDLogger gmx_unused& mdlog,
309                                           const gmx::HardwareTopology gmx_unused& hardwareTopology)
310 {
311 #if defined HAVE_SYSCONF && defined(_SC_NPROCESSORS_CONF)
312     if (hardwareTopology.supportLevel() < gmx::HardwareTopology::SupportLevel::LogicalProcessorCount)
313     {
314         return;
315     }
316
317     int countFromDetection = hardwareTopology.machine().logicalProcessorCount;
318     int countConfigured    = sysconf(_SC_NPROCESSORS_CONF);
319
320     /* BIOS, kernel or user actions can take physical processors
321      * offline. We already cater for the some of the cases inside the hardwareToplogy
322      * by trying to spin up cores just before we detect, but there could be other
323      * cases where it is worthwhile to hint that there might be more resources available.
324      */
325     if (countConfigured >= 0 && countConfigured != countFromDetection)
326     {
327         GMX_LOG(mdlog.info)
328                 .appendTextFormatted(
329                         "Note: %d CPUs configured, but only %d were detected to be online.\n",
330                         countConfigured,
331                         countFromDetection);
332
333         if (c_architecture == Architecture::X86 && countConfigured == 2 * countFromDetection)
334         {
335             GMX_LOG(mdlog.info)
336                     .appendText(
337                             "      X86 Hyperthreading is likely disabled; enable it for better "
338                             "performance.");
339         }
340         // For PowerPC (likely Power8) it is possible to set SMT to either 2,4, or 8-way hardware threads.
341         // We only warn if it is completely disabled since default performance drops with SMT8.
342         if (c_architecture == Architecture::PowerPC && countConfigured == 8 * countFromDetection)
343         {
344             GMX_LOG(mdlog.info)
345                     .appendText(
346                             "      PowerPC SMT is likely disabled; enable SMT2/SMT4 for better "
347                             "performance.");
348         }
349     }
350 #else
351     GMX_UNUSED_VALUE(mdlog);
352     GMX_UNUSED_VALUE(hardwareTopology);
353 #endif
354 }
355
356 std::unique_ptr<gmx_hw_info_t> gmx_detect_hardware(const PhysicalNodeCommunicator& physicalNodeComm)
357 {
358     // Ensure all cores have spun up, where applicable.
359     hardwareTopologyPrepareDetection();
360
361     // TODO: We should also do CPU hardware detection only once on each
362     // physical node and broadcast it, instead of doing it on every MPI rank.
363     auto hardwareInfo = std::make_unique<gmx_hw_info_t>(
364             std::make_unique<CpuInfo>(CpuInfo::detect()),
365             std::make_unique<HardwareTopology>(HardwareTopology::detect()));
366
367     // TODO: Get rid of this altogether.
368     hardwareInfo->nthreads_hw_avail = hardwareInfo->hardwareTopology->machine().logicalProcessorCount;
369
370     // Detect GPUs
371     // Open a nested scope so no temporary variables can
372     // be mis-used later.
373     {
374         DeviceDetectionResult deviceDetectionResult = detectAllDeviceInformation(physicalNodeComm);
375         hardwareInfo->deviceInfoList.swap(deviceDetectionResult.deviceInfoList_);
376         std::swap(hardwareInfo->hardwareDetectionWarnings_, deviceDetectionResult.deviceDetectionWarnings_);
377     }
378
379     gmx_collect_hardware_mpi(*hardwareInfo->cpuInfo, physicalNodeComm, hardwareInfo.get());
380
381     return hardwareInfo;
382 }
383
384 void logHardwareDetectionWarnings(const gmx::MDLogger& mdlog, const gmx_hw_info_t& hardwareInformation)
385 {
386     for (const std::string& warningString : hardwareInformation.hardwareDetectionWarnings_)
387     {
388         GMX_LOG(mdlog.warning).asParagraph().appendText(warningString);
389     }
390 }
391
392 } // namespace gmx