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