aea9a1a3e6fef03dc48374eac44285121de5b85a
[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,2017,2018, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "detecthardware.h"
38
39 #include "config.h"
40
41 #include <cerrno>
42 #include <cstdlib>
43 #include <cstring>
44
45 #include <algorithm>
46 #include <array>
47 #include <chrono>
48 #include <memory>
49 #include <string>
50 #include <thread>
51 #include <vector>
52
53 #include "thread_mpi/threads.h"
54
55 #include "gromacs/compat/make_unique.h"
56 #include "gromacs/gpu_utils/gpu_utils.h"
57 #include "gromacs/hardware/cpuinfo.h"
58 #include "gromacs/hardware/hardwaretopology.h"
59 #include "gromacs/hardware/hw_info.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/simd/support.h"
62 #include "gromacs/utility/basedefinitions.h"
63 #include "gromacs/utility/basenetwork.h"
64 #include "gromacs/utility/baseversion.h"
65 #include "gromacs/utility/cstringutil.h"
66 #include "gromacs/utility/exceptions.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/gmxassert.h"
69 #include "gromacs/utility/logger.h"
70 #include "gromacs/utility/physicalnodecommunicator.h"
71 #include "gromacs/utility/programcontext.h"
72 #include "gromacs/utility/smalloc.h"
73 #include "gromacs/utility/stringutil.h"
74 #include "gromacs/utility/sysinfo.h"
75
76 #include "architecture.h"
77
78 #ifdef HAVE_UNISTD_H
79 #    include <unistd.h>       // sysconf()
80 #endif
81
82 namespace gmx
83 {
84
85 //! Convenience macro to help us avoid ifdefs each time we use sysconf
86 #if !defined(_SC_NPROCESSORS_ONLN) && defined(_SC_NPROC_ONLN)
87 #    define _SC_NPROCESSORS_ONLN _SC_NPROC_ONLN
88 #endif
89
90 //! Convenience macro to help us avoid ifdefs each time we use sysconf
91 #if !defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROC_CONF)
92 #    define _SC_NPROCESSORS_CONF _SC_NPROC_CONF
93 #endif
94
95 //! Constant used to help minimize preprocessed code
96 static const bool bGPUBinary     = GMX_GPU != GMX_GPU_NONE;
97
98 /*! \brief The hwinfo structure (common to all threads in this process).
99  *
100  * \todo This should become a shared_ptr owned by e.g. Mdrunner::runner()
101  * that is shared across any threads as needed (e.g. for thread-MPI). That
102  * offers about the same run time performance as we get here, and avoids a
103  * lot of custom code.
104  */
105 static std::unique_ptr<gmx_hw_info_t> hwinfo_g;
106 //! A reference counter for the hwinfo structure
107 static int                            n_hwinfo = 0;
108 //! A lock to protect the hwinfo structure
109 static tMPI_Thread_mutex_t            hw_info_lock = TMPI_THREAD_MUTEX_INITIALIZER;
110
111 //! Detect GPUs, if that makes sense to attempt.
112 static void gmx_detect_gpus(const gmx::MDLogger            &mdlog,
113                             const PhysicalNodeCommunicator &physicalNodeComm)
114 {
115     hwinfo_g->gpu_info.bDetectGPUs =
116         (bGPUBinary && getenv("GMX_DISABLE_GPU_DETECTION") == nullptr);
117     if (!hwinfo_g->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 reference counting on the results of the detection to
128     // enforce it. But we can assert that this is true.
129     GMX_RELEASE_ASSERT(n_hwinfo == 0, "Cannot run GPU detection on non-master thread-MPI ranks");
130     GMX_UNUSED_VALUE(physicalNodeComm);
131     isMasterRankOfPhysicalNode = true;
132 #endif
133
134     /* The OpenCL support requires us to run detection on all ranks.
135      * With CUDA we don't need to, and prefer to detect on one rank
136      * and send the information to the other ranks over MPI. */
137     bool allRanksMustDetectGpus = (GMX_GPU == GMX_GPU_OPENCL);
138     bool gpusCanBeDetected      = false;
139     if (isMasterRankOfPhysicalNode || allRanksMustDetectGpus)
140     {
141         std::string errorMessage;
142         gpusCanBeDetected = canDetectGpus(&errorMessage);
143         if (!gpusCanBeDetected)
144         {
145             GMX_LOG(mdlog.info).asParagraph().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(&hwinfo_g->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(&hwinfo_g->gpu_info.n_dev, 1, MPI_INT, 0, physicalNodeComm.comm_);
165
166         if (hwinfo_g->gpu_info.n_dev > 0)
167         {
168             int dev_size;
169
170             dev_size = hwinfo_g->gpu_info.n_dev*sizeof_gpu_dev_info();
171
172             if (!isMasterRankOfPhysicalNode)
173             {
174                 hwinfo_g->gpu_info.gpu_dev =
175                     (struct gmx_device_info_t *)malloc(dev_size);
176             }
177             MPI_Bcast(hwinfo_g->gpu_info.gpu_dev, dev_size, MPI_BYTE,
178                       0, physicalNodeComm.comm_);
179             MPI_Bcast(&hwinfo_g->gpu_info.n_dev_compatible, 1, MPI_INT,
180                       0, physicalNodeComm.comm_);
181         }
182     }
183 #endif
184 }
185
186 //! Reduce the locally collected \p hwinfo_g over MPI ranks
187 static void gmx_collect_hardware_mpi(const gmx::CpuInfo             &cpuInfo,
188                                      const PhysicalNodeCommunicator &physicalNodeComm)
189 {
190     const int  ncore        = hwinfo_g->hardwareTopology->numberOfCores();
191     /* Zen has family=23, for now we treat future AMD CPUs like Zen */
192     const bool cpuIsAmdZen  = (cpuInfo.vendor() == CpuInfo::Vendor::Amd &&
193                                cpuInfo.family() >= 23);
194
195 #if GMX_LIB_MPI
196     int       nhwthread, ngpu, i;
197     int       gpu_hash;
198
199     nhwthread = hwinfo_g->nthreads_hw_avail;
200     ngpu      = hwinfo_g->gpu_info.n_dev_compatible;
201     /* Create a unique hash of the GPU type(s) in this node */
202     gpu_hash  = 0;
203     /* Here it might be better to only loop over the compatible GPU, but we
204      * don't have that information available and it would also require
205      * removing the device ID from the device info string.
206      */
207     for (i = 0; i < hwinfo_g->gpu_info.n_dev; i++)
208     {
209         char stmp[STRLEN];
210
211         /* Since the device ID is incorporated in the hash, the order of
212          * the GPUs affects the hash. Also two identical GPUs won't give
213          * a gpu_hash of zero after XORing.
214          */
215         get_gpu_device_info_string(stmp, hwinfo_g->gpu_info, i);
216         gpu_hash ^= gmx_string_fullhash_func(stmp, gmx_string_hash_init);
217     }
218
219     constexpr int                          numElementsCounts =  4;
220     std::array<int, numElementsCounts>     countsReduced;
221     {
222         std::array<int, numElementsCounts> countsLocal = {{0}};
223         // Organize to sum values from only one rank within each node,
224         // so we get the sum over all nodes.
225         bool isMasterRankOfPhysicalNode = (physicalNodeComm.rank_ == 0);
226         if (isMasterRankOfPhysicalNode)
227         {
228             countsLocal[0] = 1;
229             countsLocal[1] = ncore;
230             countsLocal[2] = nhwthread;
231             countsLocal[3] = ngpu;
232         }
233
234         MPI_Allreduce(countsLocal.data(), countsReduced.data(), countsLocal.size(),
235                       MPI_INT, MPI_SUM, MPI_COMM_WORLD);
236     }
237
238     constexpr int                       numElementsMax = 11;
239     std::array<int, numElementsMax>     maxMinReduced;
240     {
241         std::array<int, numElementsMax> maxMinLocal;
242         /* Store + and - values for all ranks,
243          * so we can get max+min with one MPI call.
244          */
245         maxMinLocal[0]  = ncore;
246         maxMinLocal[1]  = nhwthread;
247         maxMinLocal[2]  = ngpu;
248         maxMinLocal[3]  = static_cast<int>(gmx::simdSuggested(cpuInfo));
249         maxMinLocal[4]  = gpu_hash;
250         maxMinLocal[5]  = -maxMinLocal[0];
251         maxMinLocal[6]  = -maxMinLocal[1];
252         maxMinLocal[7]  = -maxMinLocal[2];
253         maxMinLocal[8]  = -maxMinLocal[3];
254         maxMinLocal[9]  = -maxMinLocal[4];
255         maxMinLocal[10] = (cpuIsAmdZen ? 1 : 0);
256
257         MPI_Allreduce(maxMinLocal.data(), maxMinReduced.data(), maxMinLocal.size(),
258                       MPI_INT, MPI_MAX, MPI_COMM_WORLD);
259     }
260
261     hwinfo_g->nphysicalnode       = countsReduced[0];
262     hwinfo_g->ncore_tot           = countsReduced[1];
263     hwinfo_g->ncore_min           = -maxMinReduced[5];
264     hwinfo_g->ncore_max           = maxMinReduced[0];
265     hwinfo_g->nhwthread_tot       = countsReduced[2];
266     hwinfo_g->nhwthread_min       = -maxMinReduced[6];
267     hwinfo_g->nhwthread_max       = maxMinReduced[1];
268     hwinfo_g->ngpu_compatible_tot = countsReduced[3];
269     hwinfo_g->ngpu_compatible_min = -maxMinReduced[7];
270     hwinfo_g->ngpu_compatible_max = maxMinReduced[2];
271     hwinfo_g->simd_suggest_min    = -maxMinReduced[8];
272     hwinfo_g->simd_suggest_max    = maxMinReduced[3];
273     hwinfo_g->bIdenticalGPUs      = (maxMinReduced[4] == -maxMinReduced[9]);
274     hwinfo_g->haveAmdZenCpu       = (maxMinReduced[10] > 0);
275 #else
276     /* All ranks use the same pointer, protected by a mutex in the caller */
277     hwinfo_g->nphysicalnode       = 1;
278     hwinfo_g->ncore_tot           = ncore;
279     hwinfo_g->ncore_min           = ncore;
280     hwinfo_g->ncore_max           = ncore;
281     hwinfo_g->nhwthread_tot       = hwinfo_g->nthreads_hw_avail;
282     hwinfo_g->nhwthread_min       = hwinfo_g->nthreads_hw_avail;
283     hwinfo_g->nhwthread_max       = hwinfo_g->nthreads_hw_avail;
284     hwinfo_g->ngpu_compatible_tot = hwinfo_g->gpu_info.n_dev_compatible;
285     hwinfo_g->ngpu_compatible_min = hwinfo_g->gpu_info.n_dev_compatible;
286     hwinfo_g->ngpu_compatible_max = hwinfo_g->gpu_info.n_dev_compatible;
287     hwinfo_g->simd_suggest_min    = static_cast<int>(simdSuggested(cpuInfo));
288     hwinfo_g->simd_suggest_max    = static_cast<int>(simdSuggested(cpuInfo));
289     hwinfo_g->bIdenticalGPUs      = TRUE;
290     hwinfo_g->haveAmdZenCpu       = cpuIsAmdZen;
291     GMX_UNUSED_VALUE(physicalNodeComm);
292 #endif
293 }
294
295 /*! \brief Utility that does dummy computing for max 2 seconds to spin up cores
296  *
297  *  This routine will check the number of cores configured and online
298  *  (using sysconf), and the spins doing dummy compute operations for up to
299  *  2 seconds, or until all cores have come online. This can be used prior to
300  *  hardware detection for platforms that take unused processors offline.
301  *
302  *  This routine will not throw exceptions.
303  */
304 static void
305 spinUpCore() noexcept
306 {
307 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROCESSORS_ONLN)
308     float dummy           = 0.1;
309     int   countConfigured = sysconf(_SC_NPROCESSORS_CONF);    // noexcept
310     auto  start           = std::chrono::steady_clock::now(); // noexcept
311
312     while (sysconf(_SC_NPROCESSORS_ONLN) < countConfigured &&
313            std::chrono::steady_clock::now() - start < std::chrono::seconds(2))
314     {
315         for (int i = 1; i < 10000; i++)
316         {
317             dummy /= i;
318         }
319     }
320
321     if (dummy < 0)
322     {
323         printf("This cannot happen, but prevents loop from being optimized away.");
324     }
325 #endif
326 }
327
328 /*! \brief Prepare the system before hardware topology detection
329  *
330  * This routine should perform any actions we want to put the system in a state
331  * where we want it to be before detecting the hardware topology. For most
332  * processors there is nothing to do, but some architectures (in particular ARM)
333  * have support for taking configured cores offline, which will make them disappear
334  * from the online processor count.
335  *
336  * This routine checks if there is a mismatch between the number of cores
337  * configured and online, and in that case we issue a small workload that
338  * attempts to wake sleeping cores before doing the actual detection.
339  *
340  * This type of mismatch can also occur for x86 or PowerPC on Linux, if SMT has only
341  * been disabled in the kernel (rather than bios). Since those cores will never
342  * come online automatically, we currently skip this test for x86 & PowerPC to
343  * avoid wasting 2 seconds. We also skip the test if there is no thread support.
344  *
345  * \note Cores will sleep relatively quickly again, so it's important to issue
346  *       the real detection code directly after this routine.
347  */
348 static void
349 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 &&
356         c_architecture != Architecture::PowerPC)
357     {
358         int                      countConfigured  = sysconf(_SC_NPROCESSORS_CONF);
359         std::vector<std::thread> workThreads(countConfigured);
360
361         for (auto &t : workThreads)
362         {
363             t = std::thread(spinUpCore);
364         }
365
366         for (auto &t : workThreads)
367         {
368             t.join();
369         }
370     }
371 #endif
372 }
373
374 /*! \brief Sanity check hardware topology and print some notes to log
375  *
376  *  \param mdlog            Logger.
377  *  \param hardwareTopology Reference to hardwareTopology object.
378  */
379 static void
380 hardwareTopologyDoubleCheckDetection(const gmx::MDLogger gmx_unused         &mdlog,
381                                      const gmx::HardwareTopology gmx_unused &hardwareTopology)
382 {
383 #if defined HAVE_SYSCONF && defined(_SC_NPROCESSORS_CONF)
384     if (hardwareTopology.supportLevel() < gmx::HardwareTopology::SupportLevel::LogicalProcessorCount)
385     {
386         return;
387     }
388
389     int countFromDetection = hardwareTopology.machine().logicalProcessorCount;
390     int countConfigured    = sysconf(_SC_NPROCESSORS_CONF);
391
392     /* BIOS, kernel or user actions can take physical processors
393      * offline. We already cater for the some of the cases inside the hardwareToplogy
394      * by trying to spin up cores just before we detect, but there could be other
395      * cases where it is worthwhile to hint that there might be more resources available.
396      */
397     if (countConfigured >= 0 && countConfigured != countFromDetection)
398     {
399         GMX_LOG(mdlog.info).
400             appendTextFormatted("Note: %d CPUs configured, but only %d were detected to be online.\n", countConfigured, countFromDetection);
401
402         if (c_architecture == Architecture::X86 &&
403             countConfigured == 2*countFromDetection)
404         {
405             GMX_LOG(mdlog.info).
406                 appendText("      X86 Hyperthreading is likely disabled; enable it for better performance.");
407         }
408         // For PowerPC (likely Power8) it is possible to set SMT to either 2,4, or 8-way hardware threads.
409         // We only warn if it is completely disabled since default performance drops with SMT8.
410         if (c_architecture == Architecture::PowerPC &&
411             countConfigured == 8*countFromDetection)
412         {
413             GMX_LOG(mdlog.info).
414                 appendText("      PowerPC SMT is likely disabled; enable SMT2/SMT4 for better performance.");
415         }
416     }
417 #endif
418 }
419
420 gmx_hw_info_t *gmx_detect_hardware(const gmx::MDLogger            &mdlog,
421                                    const PhysicalNodeCommunicator &physicalNodeComm)
422 {
423     int ret;
424
425     /* make sure no one else is doing the same thing */
426     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
427     if (ret != 0)
428     {
429         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
430     }
431
432     /* only initialize the hwinfo structure if it is not already initalized */
433     if (n_hwinfo == 0)
434     {
435         hwinfo_g = compat::make_unique<gmx_hw_info_t>();
436
437         /* TODO: We should also do CPU hardware detection only once on each
438          * physical node and broadcast it, instead of do it on every MPI rank. */
439         hwinfo_g->cpuInfo             = new gmx::CpuInfo(gmx::CpuInfo::detect());
440
441         hardwareTopologyPrepareDetection();
442         hwinfo_g->hardwareTopology    = new gmx::HardwareTopology(gmx::HardwareTopology::detect());
443
444         // If we detected the topology on this system, double-check that it makes sense
445         if (hwinfo_g->hardwareTopology->isThisSystem())
446         {
447             hardwareTopologyDoubleCheckDetection(mdlog, *(hwinfo_g->hardwareTopology));
448         }
449
450         // TODO: Get rid of this altogether.
451         hwinfo_g->nthreads_hw_avail = hwinfo_g->hardwareTopology->machine().logicalProcessorCount;
452
453         /* detect GPUs */
454         hwinfo_g->gpu_info.n_dev            = 0;
455         hwinfo_g->gpu_info.n_dev_compatible = 0;
456         hwinfo_g->gpu_info.gpu_dev          = nullptr;
457
458         gmx_detect_gpus(mdlog, physicalNodeComm);
459         gmx_collect_hardware_mpi(*hwinfo_g->cpuInfo, physicalNodeComm);
460     }
461     /* increase the reference counter */
462     n_hwinfo++;
463
464     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
465     if (ret != 0)
466     {
467         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
468     }
469
470     return hwinfo_g.get();
471 }
472
473 bool compatibleGpusFound(const gmx_gpu_info_t &gpu_info)
474 {
475     return gpu_info.n_dev_compatible > 0;
476 }
477
478 void gmx_hardware_info_free()
479 {
480     int ret;
481
482     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
483     if (ret != 0)
484     {
485         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
486     }
487
488     /* decrease the reference counter */
489     n_hwinfo--;
490
491
492     if (n_hwinfo < 0)
493     {
494         gmx_incons("n_hwinfo < 0");
495     }
496
497     if (n_hwinfo == 0)
498     {
499         delete hwinfo_g->cpuInfo;
500         delete hwinfo_g->hardwareTopology;
501         free_gpu_info(&hwinfo_g->gpu_info);
502         hwinfo_g.reset();
503     }
504
505     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
506     if (ret != 0)
507     {
508         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
509     }
510 }
511
512 }  // namespace gmx