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