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