Check for GPU detection support before detecting
[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, 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/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, const t_commrec *cr)
112 {
113 #if GMX_LIB_MPI
114     int              rank_world;
115     MPI_Comm         physicalnode_comm;
116 #endif
117     bool             isMasterRankOfNode;
118
119     hwinfo_g->gpu_info.bDetectGPUs =
120         (bGPUBinary && getenv("GMX_DISABLE_GPU_DETECTION") == nullptr);
121     if (!hwinfo_g->gpu_info.bDetectGPUs)
122     {
123         return;
124     }
125
126     /* Under certain circumstances MPI ranks on the same physical node
127      * can not simultaneously access the same GPU(s). Therefore we run
128      * the detection only on one MPI rank per node and broadcast the info.
129      * Note that with thread-MPI only a single thread runs this code.
130      *
131      * NOTE: We can't broadcast gpu_info with OpenCL as the device and platform
132      * ID stored in the structure are unique for each rank (even if a device
133      * is shared by multiple ranks).
134      *
135      * TODO: We should also do CPU hardware detection only once on each
136      * physical node and broadcast it, instead of do it on every MPI rank.
137      */
138 #if GMX_LIB_MPI
139     /* A split of MPI_COMM_WORLD over physical nodes is only required here,
140      * so we create and destroy it locally.
141      */
142     MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
143     MPI_Comm_split(MPI_COMM_WORLD, gmx_physicalnode_id_hash(),
144                    rank_world, &physicalnode_comm);
145     {
146         int rankOnNode = -1;
147         MPI_Comm_rank(physicalnode_comm, &rankOnNode);
148         isMasterRankOfNode = (rankOnNode == 0);
149     }
150     GMX_UNUSED_VALUE(cr);
151 #else
152     // Here there should be only one process, because if we are using
153     // thread-MPI, only one thread is active so far. So we check this.
154     GMX_RELEASE_ASSERT(cr->nnodes == 1 && cr->sim_nodeid == 0, "Only a single (master) process should execute here");
155     isMasterRankOfNode = true;
156 #endif
157
158     /*  With CUDA detect only on one rank per host, with OpenCL need do
159      *  the detection on all PP ranks */
160     bool isOpenclPpRank = ((GMX_GPU == GMX_GPU_OPENCL) && thisRankHasDuty(cr, DUTY_PP));
161
162     bool gpusCanBeDetected = false;
163     if (isMasterRankOfNode || isOpenclPpRank)
164     {
165         gpusCanBeDetected = canDetectGpus();
166         // No need to tell the user anything at this point, they get a
167         // hardware report later.
168     }
169
170     if (gpusCanBeDetected)
171     {
172         char detection_error[STRLEN] = "", sbuf[STRLEN];
173
174         if (detect_gpus(&hwinfo_g->gpu_info, detection_error) != 0)
175         {
176             if (detection_error[0] != '\0')
177             {
178                 sprintf(sbuf, ":\n      %s\n", detection_error);
179             }
180             else
181             {
182                 sprintf(sbuf, ".");
183             }
184             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
185                     "NOTE: Error occurred during GPU detection%s"
186                     "      Can not use GPU acceleration, will fall back to CPU kernels.",
187                     sbuf);
188         }
189     }
190
191 #if GMX_LIB_MPI
192     if (!isOpenclPpRank)
193     {
194         /* Broadcast the GPU info to the other ranks within this node */
195         MPI_Bcast(&hwinfo_g->gpu_info.n_dev, 1, MPI_INT, 0, physicalnode_comm);
196
197         if (hwinfo_g->gpu_info.n_dev > 0)
198         {
199             int dev_size;
200
201             dev_size = hwinfo_g->gpu_info.n_dev*sizeof_gpu_dev_info();
202
203             if (!isMasterRankOfNode)
204             {
205                 hwinfo_g->gpu_info.gpu_dev =
206                     (struct gmx_device_info_t *)malloc(dev_size);
207             }
208             MPI_Bcast(hwinfo_g->gpu_info.gpu_dev, dev_size, MPI_BYTE,
209                       0, physicalnode_comm);
210             MPI_Bcast(&hwinfo_g->gpu_info.n_dev_compatible, 1, MPI_INT,
211                       0, physicalnode_comm);
212         }
213     }
214
215     MPI_Comm_free(&physicalnode_comm);
216 #endif
217 }
218
219 //! Reduce the locally collected \p hwinfo_g over MPI ranks
220 static void gmx_collect_hardware_mpi(const gmx::CpuInfo &cpuInfo)
221 {
222     const int  ncore        = hwinfo_g->hardwareTopology->numberOfCores();
223     /* Zen has family=23, for now we treat future AMD CPUs like Zen */
224     const bool cpuIsAmdZen  = (cpuInfo.vendor() == CpuInfo::Vendor::Amd &&
225                                cpuInfo.family() >= 23);
226
227 #if GMX_LIB_MPI
228     int       rank_id;
229     int       nrank, rank, nhwthread, ngpu, i;
230     int       gpu_hash;
231     int      *buf, *all;
232
233     rank_id   = gmx_physicalnode_id_hash();
234     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
235     MPI_Comm_size(MPI_COMM_WORLD, &nrank);
236     nhwthread = hwinfo_g->nthreads_hw_avail;
237     ngpu      = hwinfo_g->gpu_info.n_dev_compatible;
238     /* Create a unique hash of the GPU type(s) in this node */
239     gpu_hash  = 0;
240     /* Here it might be better to only loop over the compatible GPU, but we
241      * don't have that information available and it would also require
242      * removing the device ID from the device info string.
243      */
244     for (i = 0; i < hwinfo_g->gpu_info.n_dev; i++)
245     {
246         char stmp[STRLEN];
247
248         /* Since the device ID is incorporated in the hash, the order of
249          * the GPUs affects the hash. Also two identical GPUs won't give
250          * a gpu_hash of zero after XORing.
251          */
252         get_gpu_device_info_string(stmp, hwinfo_g->gpu_info, i);
253         gpu_hash ^= gmx_string_fullhash_func(stmp, gmx_string_hash_init);
254     }
255
256     snew(buf, nrank);
257     snew(all, nrank);
258     buf[rank] = rank_id;
259
260     MPI_Allreduce(buf, all, nrank, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
261
262     gmx_bool bFound;
263     int      nnode0, ncore0, nhwthread0, ngpu0, r;
264
265     bFound     = FALSE;
266     ncore0     = 0;
267     nnode0     = 0;
268     nhwthread0 = 0;
269     ngpu0      = 0;
270     for (r = 0; r < nrank; r++)
271     {
272         if (all[r] == rank_id)
273         {
274             if (!bFound && r == rank)
275             {
276                 /* We are the first rank in this physical node */
277                 nnode0     = 1;
278                 ncore0     = ncore;
279                 nhwthread0 = nhwthread;
280                 ngpu0      = ngpu;
281             }
282             bFound = TRUE;
283         }
284     }
285
286     sfree(buf);
287     sfree(all);
288
289     constexpr int                          numElementsCounts =  4;
290     std::array<int, numElementsCounts>     countsReduced;
291     {
292         std::array<int, numElementsCounts> countsLocal;
293         /* Sum values from only intra-rank 0 so we get the sum over all nodes */
294         countsLocal[0] = nnode0;
295         countsLocal[1] = ncore0;
296         countsLocal[2] = nhwthread0;
297         countsLocal[3] = ngpu0;
298
299         MPI_Allreduce(countsLocal.data(), countsReduced.data(), countsLocal.size(),
300                       MPI_INT, MPI_SUM, MPI_COMM_WORLD);
301     }
302
303     constexpr int                       numElementsMax = 11;
304     std::array<int, numElementsMax>     maxMinReduced;
305     {
306         std::array<int, numElementsMax> maxMinLocal;
307         /* Store + and - values for all ranks,
308          * so we can get max+min with one MPI call.
309          */
310         maxMinLocal[0]  = ncore;
311         maxMinLocal[1]  = nhwthread;
312         maxMinLocal[2]  = ngpu;
313         maxMinLocal[3]  = static_cast<int>(gmx::simdSuggested(cpuInfo));
314         maxMinLocal[4]  = gpu_hash;
315         maxMinLocal[5]  = -maxMinLocal[0];
316         maxMinLocal[6]  = -maxMinLocal[1];
317         maxMinLocal[7]  = -maxMinLocal[2];
318         maxMinLocal[8]  = -maxMinLocal[3];
319         maxMinLocal[9]  = -maxMinLocal[4];
320         maxMinLocal[10] = (cpuIsAmdZen ? 1 : 0);
321
322         MPI_Allreduce(maxMinLocal.data(), maxMinReduced.data(), maxMinLocal.size(),
323                       MPI_INT, MPI_MAX, MPI_COMM_WORLD);
324     }
325
326     hwinfo_g->nphysicalnode       = countsReduced[0];
327     hwinfo_g->ncore_tot           = countsReduced[1];
328     hwinfo_g->ncore_min           = -maxMinReduced[5];
329     hwinfo_g->ncore_max           = maxMinReduced[0];
330     hwinfo_g->nhwthread_tot       = countsReduced[2];
331     hwinfo_g->nhwthread_min       = -maxMinReduced[6];
332     hwinfo_g->nhwthread_max       = maxMinReduced[1];
333     hwinfo_g->ngpu_compatible_tot = countsReduced[3];
334     hwinfo_g->ngpu_compatible_min = -maxMinReduced[7];
335     hwinfo_g->ngpu_compatible_max = maxMinReduced[2];
336     hwinfo_g->simd_suggest_min    = -maxMinReduced[8];
337     hwinfo_g->simd_suggest_max    = maxMinReduced[3];
338     hwinfo_g->bIdenticalGPUs      = (maxMinReduced[4] == -maxMinReduced[9]);
339     hwinfo_g->haveAmdZenCpu       = (maxMinReduced[10] > 0);
340 #else
341     /* All ranks use the same pointer, protected by a mutex in the caller */
342     hwinfo_g->nphysicalnode       = 1;
343     hwinfo_g->ncore_tot           = ncore;
344     hwinfo_g->ncore_min           = ncore;
345     hwinfo_g->ncore_max           = ncore;
346     hwinfo_g->nhwthread_tot       = hwinfo_g->nthreads_hw_avail;
347     hwinfo_g->nhwthread_min       = hwinfo_g->nthreads_hw_avail;
348     hwinfo_g->nhwthread_max       = hwinfo_g->nthreads_hw_avail;
349     hwinfo_g->ngpu_compatible_tot = hwinfo_g->gpu_info.n_dev_compatible;
350     hwinfo_g->ngpu_compatible_min = hwinfo_g->gpu_info.n_dev_compatible;
351     hwinfo_g->ngpu_compatible_max = hwinfo_g->gpu_info.n_dev_compatible;
352     hwinfo_g->simd_suggest_min    = static_cast<int>(simdSuggested(cpuInfo));
353     hwinfo_g->simd_suggest_max    = static_cast<int>(simdSuggested(cpuInfo));
354     hwinfo_g->bIdenticalGPUs      = TRUE;
355     hwinfo_g->haveAmdZenCpu       = cpuIsAmdZen;
356 #endif
357 }
358
359 /*! \brief Utility that does dummy computing for max 2 seconds to spin up cores
360  *
361  *  This routine will check the number of cores configured and online
362  *  (using sysconf), and the spins doing dummy compute operations for up to
363  *  2 seconds, or until all cores have come online. This can be used prior to
364  *  hardware detection for platforms that take unused processors offline.
365  *
366  *  This routine will not throw exceptions.
367  */
368 static void
369 spinUpCore() noexcept
370 {
371 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROCESSORS_ONLN)
372     float dummy           = 0.1;
373     int   countConfigured = sysconf(_SC_NPROCESSORS_CONF);    // noexcept
374     auto  start           = std::chrono::steady_clock::now(); // noexcept
375
376     while (sysconf(_SC_NPROCESSORS_ONLN) < countConfigured &&
377            std::chrono::steady_clock::now() - start < std::chrono::seconds(2))
378     {
379         for (int i = 1; i < 10000; i++)
380         {
381             dummy /= i;
382         }
383     }
384
385     if (dummy < 0)
386     {
387         printf("This cannot happen, but prevents loop from being optimized away.");
388     }
389 #endif
390 }
391
392 /*! \brief Prepare the system before hardware topology detection
393  *
394  * This routine should perform any actions we want to put the system in a state
395  * where we want it to be before detecting the hardware topology. For most
396  * processors there is nothing to do, but some architectures (in particular ARM)
397  * have support for taking configured cores offline, which will make them disappear
398  * from the online processor count.
399  *
400  * This routine checks if there is a mismatch between the number of cores
401  * configured and online, and in that case we issue a small workload that
402  * attempts to wake sleeping cores before doing the actual detection.
403  *
404  * This type of mismatch can also occur for x86 or PowerPC on Linux, if SMT has only
405  * been disabled in the kernel (rather than bios). Since those cores will never
406  * come online automatically, we currently skip this test for x86 & PowerPC to
407  * avoid wasting 2 seconds. We also skip the test if there is no thread support.
408  *
409  * \note Cores will sleep relatively quickly again, so it's important to issue
410  *       the real detection code directly after this routine.
411  */
412 static void
413 hardwareTopologyPrepareDetection()
414 {
415 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && \
416     (defined(THREAD_PTHREADS) || defined(THREAD_WINDOWS))
417
418     // Modify this conditional when/if x86 or PowerPC starts to sleep some cores
419     if (c_architecture != Architecture::X86 &&
420         c_architecture != Architecture::PowerPC)
421     {
422         int                      countConfigured  = sysconf(_SC_NPROCESSORS_CONF);
423         std::vector<std::thread> workThreads(countConfigured);
424
425         for (auto &t : workThreads)
426         {
427             t = std::thread(spinUpCore);
428         }
429
430         for (auto &t : workThreads)
431         {
432             t.join();
433         }
434     }
435 #endif
436 }
437
438 /*! \brief Sanity check hardware topology and print some notes to log
439  *
440  *  \param mdlog            Logger.
441  *  \param hardwareTopology Reference to hardwareTopology object.
442  */
443 static void
444 hardwareTopologyDoubleCheckDetection(const gmx::MDLogger gmx_unused         &mdlog,
445                                      const gmx::HardwareTopology gmx_unused &hardwareTopology)
446 {
447 #if defined HAVE_SYSCONF && defined(_SC_NPROCESSORS_CONF)
448     if (hardwareTopology.supportLevel() < gmx::HardwareTopology::SupportLevel::LogicalProcessorCount)
449     {
450         return;
451     }
452
453     int countFromDetection = hardwareTopology.machine().logicalProcessorCount;
454     int countConfigured    = sysconf(_SC_NPROCESSORS_CONF);
455
456     /* BIOS, kernel or user actions can take physical processors
457      * offline. We already cater for the some of the cases inside the hardwareToplogy
458      * by trying to spin up cores just before we detect, but there could be other
459      * cases where it is worthwhile to hint that there might be more resources available.
460      */
461     if (countConfigured >= 0 && countConfigured != countFromDetection)
462     {
463         GMX_LOG(mdlog.info).
464             appendTextFormatted("Note: %d CPUs configured, but only %d were detected to be online.\n", countConfigured, countFromDetection);
465
466         if (c_architecture == Architecture::X86 &&
467             countConfigured == 2*countFromDetection)
468         {
469             GMX_LOG(mdlog.info).
470                 appendText("      X86 Hyperthreading is likely disabled; enable it for better performance.");
471         }
472         // For PowerPC (likely Power8) it is possible to set SMT to either 2,4, or 8-way hardware threads.
473         // We only warn if it is completely disabled since default performance drops with SMT8.
474         if (c_architecture == Architecture::PowerPC &&
475             countConfigured == 8*countFromDetection)
476         {
477             GMX_LOG(mdlog.info).
478                 appendText("      PowerPC SMT is likely disabled; enable SMT2/SMT4 for better performance.");
479         }
480     }
481 #endif
482 }
483
484 gmx_hw_info_t *gmx_detect_hardware(const gmx::MDLogger &mdlog, const t_commrec *cr)
485 {
486     int ret;
487
488     /* make sure no one else is doing the same thing */
489     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
490     if (ret != 0)
491     {
492         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
493     }
494
495     /* only initialize the hwinfo structure if it is not already initalized */
496     if (n_hwinfo == 0)
497     {
498         hwinfo_g = compat::make_unique<gmx_hw_info_t>();
499
500         hwinfo_g->cpuInfo             = new gmx::CpuInfo(gmx::CpuInfo::detect());
501
502         hardwareTopologyPrepareDetection();
503         hwinfo_g->hardwareTopology    = new gmx::HardwareTopology(gmx::HardwareTopology::detect());
504
505         // If we detected the topology on this system, double-check that it makes sense
506         if (hwinfo_g->hardwareTopology->isThisSystem())
507         {
508             hardwareTopologyDoubleCheckDetection(mdlog, *(hwinfo_g->hardwareTopology));
509         }
510
511         // TODO: Get rid of this altogether.
512         hwinfo_g->nthreads_hw_avail = hwinfo_g->hardwareTopology->machine().logicalProcessorCount;
513
514         /* detect GPUs */
515         hwinfo_g->gpu_info.n_dev            = 0;
516         hwinfo_g->gpu_info.n_dev_compatible = 0;
517         hwinfo_g->gpu_info.gpu_dev          = nullptr;
518
519         gmx_detect_gpus(mdlog, cr);
520         gmx_collect_hardware_mpi(*hwinfo_g->cpuInfo);
521     }
522     /* increase the reference counter */
523     n_hwinfo++;
524
525     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
526     if (ret != 0)
527     {
528         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
529     }
530
531     return hwinfo_g.get();
532 }
533
534 bool compatibleGpusFound(const gmx_gpu_info_t &gpu_info)
535 {
536     return gpu_info.n_dev_compatible > 0;
537 }
538
539 void gmx_hardware_info_free()
540 {
541     int ret;
542
543     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
544     if (ret != 0)
545     {
546         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
547     }
548
549     /* decrease the reference counter */
550     n_hwinfo--;
551
552
553     if (n_hwinfo < 0)
554     {
555         gmx_incons("n_hwinfo < 0");
556     }
557
558     if (n_hwinfo == 0)
559     {
560         delete hwinfo_g->cpuInfo;
561         delete hwinfo_g->hardwareTopology;
562         free_gpu_info(&hwinfo_g->gpu_info);
563         hwinfo_g.reset();
564     }
565
566     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
567     if (ret != 0)
568     {
569         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
570     }
571 }
572
573 }  // namespace gmx