Merge remote-tracking branch 'origin/release-2019'
[gromacs.git] / src / gromacs / hardware / detecthardware.cpp
blob41071a68725c9f169ac53c72d74ebe818544a6b5
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, 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.
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.
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.
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.
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.
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.
35 #include "gmxpre.h"
37 #include "detecthardware.h"
39 #include "config.h"
41 #include <algorithm>
42 #include <array>
43 #include <chrono>
44 #include <memory>
45 #include <string>
46 #include <thread>
47 #include <vector>
49 #include "gromacs/compat/pointers.h"
50 #include "gromacs/gpu_utils/gpu_utils.h"
51 #include "gromacs/hardware/cpuinfo.h"
52 #include "gromacs/hardware/hardwaretopology.h"
53 #include "gromacs/hardware/hw_info.h"
54 #include "gromacs/simd/support.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/basenetwork.h"
57 #include "gromacs/utility/baseversion.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/exceptions.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/logger.h"
63 #include "gromacs/utility/mutex.h"
64 #include "gromacs/utility/physicalnodecommunicator.h"
66 #include "architecture.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()
81 free_gpu_info(&gpu_info);
84 namespace gmx
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
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
97 /*! \brief Information about the hardware of all nodes (common to all threads in this process).
99 * This information is constructed only when required, but thereafter
100 * its lifetime is that of the whole process, potentially across
101 * multiple successive simulation parts. It's wise to ensure that only
102 * one thread can create the information, but thereafter they can all
103 * read it without e.g. needing a std::shared_ptr to ensure its
104 * lifetime exceeds that of the thread. */
105 static std::unique_ptr<gmx_hw_info_t> g_hardwareInfo;
106 //! A mutex to protect the hwinfo structure
107 static Mutex g_hardwareInfoMutex;
109 //! Detect GPUs, if that makes sense to attempt.
110 static void gmx_detect_gpus(const gmx::MDLogger &mdlog,
111 const PhysicalNodeCommunicator &physicalNodeComm,
112 compat::not_null<gmx_hw_info_t *> hardwareInfo)
114 hardwareInfo->gpu_info.bDetectGPUs = canPerformGpuDetection();
116 if (!hardwareInfo->gpu_info.bDetectGPUs)
118 return;
121 bool isMasterRankOfPhysicalNode = true;
122 #if GMX_LIB_MPI
123 isMasterRankOfPhysicalNode = (physicalNodeComm.rank_ == 0);
124 #else
125 // We choose to run the detection only once with thread-MPI and
126 // use a mutex to enforce it.
127 GMX_UNUSED_VALUE(physicalNodeComm);
128 isMasterRankOfPhysicalNode = true;
129 #endif
131 /* The OpenCL support requires us to run detection on all ranks.
132 * With CUDA we don't need to, and prefer to detect on one rank
133 * and send the information to the other ranks over MPI. */
134 bool allRanksMustDetectGpus = (GMX_GPU == GMX_GPU_OPENCL);
135 bool gpusCanBeDetected = false;
136 if (isMasterRankOfPhysicalNode || allRanksMustDetectGpus)
138 std::string errorMessage;
139 gpusCanBeDetected = isGpuDetectionFunctional(&errorMessage);
140 if (!gpusCanBeDetected)
142 GMX_LOG(mdlog.info).asParagraph().appendTextFormatted(
143 "NOTE: Detection of GPUs failed. The API reported:\n"
144 " %s\n"
145 " GROMACS cannot run tasks on a GPU.",
146 errorMessage.c_str());
150 if (gpusCanBeDetected)
152 findGpus(&hardwareInfo->gpu_info);
153 // No need to tell the user anything at this point, they get a
154 // hardware report later.
157 #if GMX_LIB_MPI
158 if (!allRanksMustDetectGpus)
160 /* Broadcast the GPU info to the other ranks within this node */
161 MPI_Bcast(&hardwareInfo->gpu_info.n_dev, 1, MPI_INT, 0, physicalNodeComm.comm_);
163 if (hardwareInfo->gpu_info.n_dev > 0)
165 int dev_size;
167 dev_size = hardwareInfo->gpu_info.n_dev*sizeof_gpu_dev_info();
169 if (!isMasterRankOfPhysicalNode)
171 hardwareInfo->gpu_info.gpu_dev =
172 (struct gmx_device_info_t *)malloc(dev_size);
174 MPI_Bcast(hardwareInfo->gpu_info.gpu_dev, dev_size, MPI_BYTE,
175 0, physicalNodeComm.comm_);
176 MPI_Bcast(&hardwareInfo->gpu_info.n_dev_compatible, 1, MPI_INT,
177 0, physicalNodeComm.comm_);
180 #endif
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)
188 const int ncore = hardwareInfo->hardwareTopology->numberOfCores();
189 /* Zen1 is assumed for:
190 * - family=23 with the below listed models;
191 * - Hygon as vendor.
193 const bool cpuIsAmdZen1 = ((cpuInfo.vendor() == CpuInfo::Vendor::Amd &&
194 cpuInfo.family() == 23 &&
195 (cpuInfo.model() == 1 || cpuInfo.model() == 17 ||
196 cpuInfo.model() == 8 || cpuInfo.model() == 24)) ||
197 cpuInfo.vendor() == CpuInfo::Vendor::Hygon);
198 #if GMX_LIB_MPI
199 int nhwthread, ngpu, i;
200 int gpu_hash;
202 nhwthread = hardwareInfo->nthreads_hw_avail;
203 ngpu = hardwareInfo->gpu_info.n_dev_compatible;
204 /* Create a unique hash of the GPU type(s) in this node */
205 gpu_hash = 0;
206 /* Here it might be better to only loop over the compatible GPU, but we
207 * don't have that information available and it would also require
208 * removing the device ID from the device info string.
210 for (i = 0; i < hardwareInfo->gpu_info.n_dev; i++)
212 char stmp[STRLEN];
214 /* Since the device ID is incorporated in the hash, the order of
215 * the GPUs affects the hash. Also two identical GPUs won't give
216 * a gpu_hash of zero after XORing.
218 get_gpu_device_info_string(stmp, hardwareInfo->gpu_info, i);
219 gpu_hash ^= gmx_string_fullhash_func(stmp, gmx_string_hash_init);
222 constexpr int numElementsCounts = 4;
223 std::array<int, numElementsCounts> countsReduced;
225 std::array<int, numElementsCounts> countsLocal = {{0}};
226 // Organize to sum values from only one rank within each node,
227 // so we get the sum over all nodes.
228 bool isMasterRankOfPhysicalNode = (physicalNodeComm.rank_ == 0);
229 if (isMasterRankOfPhysicalNode)
231 countsLocal[0] = 1;
232 countsLocal[1] = ncore;
233 countsLocal[2] = nhwthread;
234 countsLocal[3] = ngpu;
237 MPI_Allreduce(countsLocal.data(), countsReduced.data(), countsLocal.size(),
238 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
241 constexpr int numElementsMax = 11;
242 std::array<int, numElementsMax> maxMinReduced;
244 std::array<int, numElementsMax> maxMinLocal;
245 /* Store + and - values for all ranks,
246 * so we can get max+min with one MPI call.
248 maxMinLocal[0] = ncore;
249 maxMinLocal[1] = nhwthread;
250 maxMinLocal[2] = ngpu;
251 maxMinLocal[3] = static_cast<int>(gmx::simdSuggested(cpuInfo));
252 maxMinLocal[4] = gpu_hash;
253 maxMinLocal[5] = -maxMinLocal[0];
254 maxMinLocal[6] = -maxMinLocal[1];
255 maxMinLocal[7] = -maxMinLocal[2];
256 maxMinLocal[8] = -maxMinLocal[3];
257 maxMinLocal[9] = -maxMinLocal[4];
258 maxMinLocal[10] = (cpuIsAmdZen1 ? 1 : 0);
260 MPI_Allreduce(maxMinLocal.data(), maxMinReduced.data(), maxMinLocal.size(),
261 MPI_INT, MPI_MAX, MPI_COMM_WORLD);
264 hardwareInfo->nphysicalnode = countsReduced[0];
265 hardwareInfo->ncore_tot = countsReduced[1];
266 hardwareInfo->ncore_min = -maxMinReduced[5];
267 hardwareInfo->ncore_max = maxMinReduced[0];
268 hardwareInfo->nhwthread_tot = countsReduced[2];
269 hardwareInfo->nhwthread_min = -maxMinReduced[6];
270 hardwareInfo->nhwthread_max = maxMinReduced[1];
271 hardwareInfo->ngpu_compatible_tot = countsReduced[3];
272 hardwareInfo->ngpu_compatible_min = -maxMinReduced[7];
273 hardwareInfo->ngpu_compatible_max = maxMinReduced[2];
274 hardwareInfo->simd_suggest_min = -maxMinReduced[8];
275 hardwareInfo->simd_suggest_max = maxMinReduced[3];
276 hardwareInfo->bIdenticalGPUs = (maxMinReduced[4] == -maxMinReduced[9]);
277 hardwareInfo->haveAmdZen1Cpu = (maxMinReduced[10] > 0);
278 #else
279 /* All ranks use the same pointer, protected by a mutex in the caller */
280 hardwareInfo->nphysicalnode = 1;
281 hardwareInfo->ncore_tot = ncore;
282 hardwareInfo->ncore_min = ncore;
283 hardwareInfo->ncore_max = ncore;
284 hardwareInfo->nhwthread_tot = hardwareInfo->nthreads_hw_avail;
285 hardwareInfo->nhwthread_min = hardwareInfo->nthreads_hw_avail;
286 hardwareInfo->nhwthread_max = hardwareInfo->nthreads_hw_avail;
287 hardwareInfo->ngpu_compatible_tot = hardwareInfo->gpu_info.n_dev_compatible;
288 hardwareInfo->ngpu_compatible_min = hardwareInfo->gpu_info.n_dev_compatible;
289 hardwareInfo->ngpu_compatible_max = hardwareInfo->gpu_info.n_dev_compatible;
290 hardwareInfo->simd_suggest_min = static_cast<int>(simdSuggested(cpuInfo));
291 hardwareInfo->simd_suggest_max = static_cast<int>(simdSuggested(cpuInfo));
292 hardwareInfo->bIdenticalGPUs = TRUE;
293 hardwareInfo->haveAmdZen1Cpu = cpuIsAmdZen1;
294 GMX_UNUSED_VALUE(physicalNodeComm);
295 #endif
298 /*! \brief Utility that does dummy computing for max 2 seconds to spin up cores
300 * This routine will check the number of cores configured and online
301 * (using sysconf), and the spins doing dummy compute operations for up to
302 * 2 seconds, or until all cores have come online. This can be used prior to
303 * hardware detection for platforms that take unused processors offline.
305 * This routine will not throw exceptions.
307 static void
308 spinUpCore() noexcept
310 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROCESSORS_ONLN)
311 float dummy = 0.1;
312 int countConfigured = sysconf(_SC_NPROCESSORS_CONF); // noexcept
313 auto start = std::chrono::steady_clock::now(); // noexcept
315 while (sysconf(_SC_NPROCESSORS_ONLN) < countConfigured &&
316 std::chrono::steady_clock::now() - start < std::chrono::seconds(2))
318 for (int i = 1; i < 10000; i++)
320 dummy /= i;
324 if (dummy < 0)
326 printf("This cannot happen, but prevents loop from being optimized away.");
328 #endif
331 /*! \brief Prepare the system before hardware topology detection
333 * This routine should perform any actions we want to put the system in a state
334 * where we want it to be before detecting the hardware topology. For most
335 * processors there is nothing to do, but some architectures (in particular ARM)
336 * have support for taking configured cores offline, which will make them disappear
337 * from the online processor count.
339 * This routine checks if there is a mismatch between the number of cores
340 * configured and online, and in that case we issue a small workload that
341 * attempts to wake sleeping cores before doing the actual detection.
343 * This type of mismatch can also occur for x86 or PowerPC on Linux, if SMT has only
344 * been disabled in the kernel (rather than bios). Since those cores will never
345 * come online automatically, we currently skip this test for x86 & PowerPC to
346 * avoid wasting 2 seconds. We also skip the test if there is no thread support.
348 * \note Cores will sleep relatively quickly again, so it's important to issue
349 * the real detection code directly after this routine.
351 static void
352 hardwareTopologyPrepareDetection()
354 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && \
355 (defined(THREAD_PTHREADS) || defined(THREAD_WINDOWS))
357 // Modify this conditional when/if x86 or PowerPC starts to sleep some cores
358 if (c_architecture != Architecture::X86 &&
359 c_architecture != Architecture::PowerPC)
361 int countConfigured = sysconf(_SC_NPROCESSORS_CONF);
362 std::vector<std::thread> workThreads(countConfigured);
364 for (auto &t : workThreads)
366 t = std::thread(spinUpCore);
369 for (auto &t : workThreads)
371 t.join();
374 #endif
377 /*! \brief Sanity check hardware topology and print some notes to log
379 * \param mdlog Logger.
380 * \param hardwareTopology Reference to hardwareTopology object.
382 static void
383 hardwareTopologyDoubleCheckDetection(const gmx::MDLogger gmx_unused &mdlog,
384 const gmx::HardwareTopology gmx_unused &hardwareTopology)
386 #if defined HAVE_SYSCONF && defined(_SC_NPROCESSORS_CONF)
387 if (hardwareTopology.supportLevel() < gmx::HardwareTopology::SupportLevel::LogicalProcessorCount)
389 return;
392 int countFromDetection = hardwareTopology.machine().logicalProcessorCount;
393 int countConfigured = sysconf(_SC_NPROCESSORS_CONF);
395 /* BIOS, kernel or user actions can take physical processors
396 * offline. We already cater for the some of the cases inside the hardwareToplogy
397 * by trying to spin up cores just before we detect, but there could be other
398 * cases where it is worthwhile to hint that there might be more resources available.
400 if (countConfigured >= 0 && countConfigured != countFromDetection)
402 GMX_LOG(mdlog.info).
403 appendTextFormatted("Note: %d CPUs configured, but only %d were detected to be online.\n", countConfigured, countFromDetection);
405 if (c_architecture == Architecture::X86 &&
406 countConfigured == 2*countFromDetection)
408 GMX_LOG(mdlog.info).
409 appendText(" X86 Hyperthreading is likely disabled; enable it for better performance.");
411 // For PowerPC (likely Power8) it is possible to set SMT to either 2,4, or 8-way hardware threads.
412 // We only warn if it is completely disabled since default performance drops with SMT8.
413 if (c_architecture == Architecture::PowerPC &&
414 countConfigured == 8*countFromDetection)
416 GMX_LOG(mdlog.info).
417 appendText(" PowerPC SMT is likely disabled; enable SMT2/SMT4 for better performance.");
420 #endif
423 gmx_hw_info_t *gmx_detect_hardware(const gmx::MDLogger &mdlog,
424 const PhysicalNodeCommunicator &physicalNodeComm)
426 // By construction, only one thread ever runs hardware detection,
427 // but we may as well prevent issues arising if that would change.
428 // Taking the lock early ensures that exactly one thread can
429 // attempt to construct g_hardwareInfo.
430 lock_guard<Mutex> lock(g_hardwareInfoMutex);
432 // If we already have the information, just return a handle to it.
433 if (g_hardwareInfo != nullptr)
435 return g_hardwareInfo.get();
438 // Make the new hardwareInfo in a temporary.
439 hardwareTopologyPrepareDetection();
441 // TODO: We should also do CPU hardware detection only once on each
442 // physical node and broadcast it, instead of doing it on every MPI rank.
443 auto hardwareInfo = std::make_unique<gmx_hw_info_t>(std::make_unique<CpuInfo>(CpuInfo::detect()),
444 std::make_unique<HardwareTopology>(HardwareTopology::detect()));
446 // If we detected the topology on this system, double-check that it makes sense
447 if (hardwareInfo->hardwareTopology->isThisSystem())
449 hardwareTopologyDoubleCheckDetection(mdlog, *hardwareInfo->hardwareTopology);
452 // TODO: Get rid of this altogether.
453 hardwareInfo->nthreads_hw_avail = hardwareInfo->hardwareTopology->machine().logicalProcessorCount;
455 // Detect GPUs
456 hardwareInfo->gpu_info.n_dev = 0;
457 hardwareInfo->gpu_info.n_dev_compatible = 0;
458 hardwareInfo->gpu_info.gpu_dev = nullptr;
460 gmx_detect_gpus(mdlog, physicalNodeComm, compat::make_not_null(hardwareInfo));
461 gmx_collect_hardware_mpi(*hardwareInfo->cpuInfo, physicalNodeComm, compat::make_not_null(hardwareInfo));
463 // Now that the temporary is fully constructed, swap it to become
464 // the real thing.
465 g_hardwareInfo.swap(hardwareInfo);
467 return g_hardwareInfo.get();
470 bool compatibleGpusFound(const gmx_gpu_info_t &gpu_info)
472 return gpu_info.n_dev_compatible > 0;
475 } // namespace gmx