Replace consistency checking with assertion
[gromacs.git] / src / gromacs / hardware / detecthardware.cpp
blob9c5bcdafd159af68124a2c4d67cdf8814e0bd31b
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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 <cerrno>
42 #include <cstdlib>
43 #include <cstring>
45 #include <algorithm>
46 #include <chrono>
47 #include <string>
48 #include <thread>
49 #include <vector>
51 #include "thread_mpi/threads.h"
53 #include "gromacs/gpu_utils/gpu_utils.h"
54 #include "gromacs/hardware/cpuinfo.h"
55 #include "gromacs/hardware/gpu_hw_info.h"
56 #include "gromacs/hardware/hardwaretopology.h"
57 #include "gromacs/hardware/hw_info.h"
58 #include "gromacs/mdtypes/commrec.h"
59 #include "gromacs/simd/support.h"
60 #include "gromacs/utility/basedefinitions.h"
61 #include "gromacs/utility/basenetwork.h"
62 #include "gromacs/utility/baseversion.h"
63 #include "gromacs/utility/cstringutil.h"
64 #include "gromacs/utility/exceptions.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/gmxassert.h"
67 #include "gromacs/utility/logger.h"
68 #include "gromacs/utility/programcontext.h"
69 #include "gromacs/utility/smalloc.h"
70 #include "gromacs/utility/stringutil.h"
71 #include "gromacs/utility/sysinfo.h"
73 #ifdef HAVE_UNISTD_H
74 # include <unistd.h> // sysconf()
75 #endif
77 //! Convenience macro to help us avoid ifdefs each time we use sysconf
78 #if !defined(_SC_NPROCESSORS_ONLN) && defined(_SC_NPROC_ONLN)
79 # define _SC_NPROCESSORS_ONLN _SC_NPROC_ONLN
80 #endif
82 //! Convenience macro to help us avoid ifdefs each time we use sysconf
83 #if !defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROC_CONF)
84 # define _SC_NPROCESSORS_CONF _SC_NPROC_CONF
85 #endif
87 #if defined (__i386__) || defined (__x86_64__) || defined (_M_IX86) || defined (_M_X64)
88 //! Constant used to help minimize preprocessed code
89 static const bool isX86 = true;
90 #else
91 //! Constant used to help minimize preprocessed code
92 static const bool isX86 = false;
93 #endif
95 #if defined __powerpc__ || defined __ppc__ || defined __PPC__
96 static const bool isPowerPC = true;
97 #else
98 static const bool isPowerPC = false;
99 #endif
101 //! Constant used to help minimize preprocessed code
102 static const bool bGPUBinary = GMX_GPU != GMX_GPU_NONE;
104 /* Note that some of the following arrays must match the "GPU support
105 * enumeration" in src/config.h.cmakein, so that GMX_GPU looks up an
106 * array entry. */
108 // TODO If/when we unify CUDA and OpenCL support code, this should
109 // move to a single place in gpu_utils.
110 /* Names of the GPU detection/check results (see e_gpu_detect_res_t in hw_info.h). */
111 const char * const gpu_detect_res_str[egpuNR] =
113 "compatible", "inexistent", "incompatible", "insane"
116 /* The globally shared hwinfo structure. */
117 static gmx_hw_info_t *hwinfo_g;
118 /* A reference counter for the hwinfo structure */
119 static int n_hwinfo = 0;
120 /* A lock to protect the hwinfo structure */
121 static tMPI_Thread_mutex_t hw_info_lock = TMPI_THREAD_MUTEX_INITIALIZER;
123 static void gmx_detect_gpus(const gmx::MDLogger &mdlog, const t_commrec *cr)
125 #if GMX_LIB_MPI
126 int rank_world;
127 MPI_Comm physicalnode_comm;
128 #endif
129 int rank_local;
131 /* Under certain circumstances MPI ranks on the same physical node
132 * can not simultaneously access the same GPU(s). Therefore we run
133 * the detection only on one MPI rank per node and broadcast the info.
134 * Note that with thread-MPI only a single thread runs this code.
136 * NOTE: We can't broadcast gpu_info with OpenCL as the device and platform
137 * ID stored in the structure are unique for each rank (even if a device
138 * is shared by multiple ranks).
140 * TODO: We should also do CPU hardware detection only once on each
141 * physical node and broadcast it, instead of do it on every MPI rank.
143 #if GMX_LIB_MPI
144 /* A split of MPI_COMM_WORLD over physical nodes is only required here,
145 * so we create and destroy it locally.
147 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
148 MPI_Comm_split(MPI_COMM_WORLD, gmx_physicalnode_id_hash(),
149 rank_world, &physicalnode_comm);
150 MPI_Comm_rank(physicalnode_comm, &rank_local);
151 GMX_UNUSED_VALUE(cr);
152 #else
153 /* Here there should be only one process, check this */
154 GMX_RELEASE_ASSERT(cr->nnodes == 1 && cr->sim_nodeid == 0, "Only a single (master) process should execute here");
156 rank_local = 0;
157 #endif
159 /* With CUDA detect only on one rank per host, with OpenCL need do
160 * the detection on all PP ranks */
161 bool isOpenclPpRank = ((GMX_GPU == GMX_GPU_OPENCL) && (cr->duty & DUTY_PP));
163 if (rank_local == 0 || isOpenclPpRank)
165 char detection_error[STRLEN] = "", sbuf[STRLEN];
167 if (detect_gpus(&hwinfo_g->gpu_info, detection_error) != 0)
169 if (detection_error[0] != '\0')
171 sprintf(sbuf, ":\n %s\n", detection_error);
173 else
175 sprintf(sbuf, ".");
177 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
178 "NOTE: Error occurred during GPU detection%s"
179 " Can not use GPU acceleration, will fall back to CPU kernels.",
180 sbuf);
184 #if GMX_LIB_MPI
185 if (!isOpenclPpRank)
187 /* Broadcast the GPU info to the other ranks within this node */
188 MPI_Bcast(&hwinfo_g->gpu_info.n_dev, 1, MPI_INT, 0, physicalnode_comm);
190 if (hwinfo_g->gpu_info.n_dev > 0)
192 int dev_size;
194 dev_size = hwinfo_g->gpu_info.n_dev*sizeof_gpu_dev_info();
196 if (rank_local > 0)
198 hwinfo_g->gpu_info.gpu_dev =
199 (struct gmx_device_info_t *)malloc(dev_size);
201 MPI_Bcast(hwinfo_g->gpu_info.gpu_dev, dev_size, MPI_BYTE,
202 0, physicalnode_comm);
203 MPI_Bcast(&hwinfo_g->gpu_info.n_dev_compatible, 1, MPI_INT,
204 0, physicalnode_comm);
208 MPI_Comm_free(&physicalnode_comm);
209 #endif
212 static void gmx_collect_hardware_mpi(const gmx::CpuInfo &cpuInfo)
214 const int ncore = hwinfo_g->hardwareTopology->numberOfCores();
215 #if GMX_LIB_MPI
216 int rank_id;
217 int nrank, rank, nhwthread, ngpu, i;
218 int gpu_hash;
219 int *buf, *all;
221 rank_id = gmx_physicalnode_id_hash();
222 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
223 MPI_Comm_size(MPI_COMM_WORLD, &nrank);
224 nhwthread = hwinfo_g->nthreads_hw_avail;
225 ngpu = hwinfo_g->gpu_info.n_dev_compatible;
226 /* Create a unique hash of the GPU type(s) in this node */
227 gpu_hash = 0;
228 /* Here it might be better to only loop over the compatible GPU, but we
229 * don't have that information available and it would also require
230 * removing the device ID from the device info string.
232 for (i = 0; i < hwinfo_g->gpu_info.n_dev; i++)
234 char stmp[STRLEN];
236 /* Since the device ID is incorporated in the hash, the order of
237 * the GPUs affects the hash. Also two identical GPUs won't give
238 * a gpu_hash of zero after XORing.
240 get_gpu_device_info_string(stmp, hwinfo_g->gpu_info, i);
241 gpu_hash ^= gmx_string_fullhash_func(stmp, gmx_string_hash_init);
244 snew(buf, nrank);
245 snew(all, nrank);
246 buf[rank] = rank_id;
248 MPI_Allreduce(buf, all, nrank, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
250 gmx_bool bFound;
251 int nnode0, ncore0, nhwthread0, ngpu0, r;
253 bFound = FALSE;
254 ncore0 = 0;
255 nnode0 = 0;
256 nhwthread0 = 0;
257 ngpu0 = 0;
258 for (r = 0; r < nrank; r++)
260 if (all[r] == rank_id)
262 if (!bFound && r == rank)
264 /* We are the first rank in this physical node */
265 nnode0 = 1;
266 ncore0 = ncore;
267 nhwthread0 = nhwthread;
268 ngpu0 = ngpu;
270 bFound = TRUE;
274 sfree(buf);
275 sfree(all);
277 int sum[4], maxmin[10];
280 int buf[4];
282 /* Sum values from only intra-rank 0 so we get the sum over all nodes */
283 buf[0] = nnode0;
284 buf[1] = ncore0;
285 buf[2] = nhwthread0;
286 buf[3] = ngpu0;
288 MPI_Allreduce(buf, sum, 4, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
292 int buf[10];
294 /* Store + and - values for all ranks,
295 * so we can get max+min with one MPI call.
297 buf[0] = ncore;
298 buf[1] = nhwthread;
299 buf[2] = ngpu;
300 buf[3] = static_cast<int>(gmx::simdSuggested(cpuInfo));
301 buf[4] = gpu_hash;
302 buf[5] = -buf[0];
303 buf[6] = -buf[1];
304 buf[7] = -buf[2];
305 buf[8] = -buf[3];
306 buf[9] = -buf[4];
308 MPI_Allreduce(buf, maxmin, 10, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
311 hwinfo_g->nphysicalnode = sum[0];
312 hwinfo_g->ncore_tot = sum[1];
313 hwinfo_g->ncore_min = -maxmin[5];
314 hwinfo_g->ncore_max = maxmin[0];
315 hwinfo_g->nhwthread_tot = sum[2];
316 hwinfo_g->nhwthread_min = -maxmin[6];
317 hwinfo_g->nhwthread_max = maxmin[1];
318 hwinfo_g->ngpu_compatible_tot = sum[3];
319 hwinfo_g->ngpu_compatible_min = -maxmin[7];
320 hwinfo_g->ngpu_compatible_max = maxmin[2];
321 hwinfo_g->simd_suggest_min = -maxmin[8];
322 hwinfo_g->simd_suggest_max = maxmin[3];
323 hwinfo_g->bIdenticalGPUs = (maxmin[4] == -maxmin[9]);
324 #else
325 /* All ranks use the same pointer, protected by a mutex in the caller */
326 hwinfo_g->nphysicalnode = 1;
327 hwinfo_g->ncore_tot = ncore;
328 hwinfo_g->ncore_min = ncore;
329 hwinfo_g->ncore_max = ncore;
330 hwinfo_g->nhwthread_tot = hwinfo_g->nthreads_hw_avail;
331 hwinfo_g->nhwthread_min = hwinfo_g->nthreads_hw_avail;
332 hwinfo_g->nhwthread_max = hwinfo_g->nthreads_hw_avail;
333 hwinfo_g->ngpu_compatible_tot = hwinfo_g->gpu_info.n_dev_compatible;
334 hwinfo_g->ngpu_compatible_min = hwinfo_g->gpu_info.n_dev_compatible;
335 hwinfo_g->ngpu_compatible_max = hwinfo_g->gpu_info.n_dev_compatible;
336 hwinfo_g->simd_suggest_min = static_cast<int>(simdSuggested(cpuInfo));
337 hwinfo_g->simd_suggest_max = static_cast<int>(simdSuggested(cpuInfo));
338 hwinfo_g->bIdenticalGPUs = TRUE;
339 #endif
342 /*! \brief Utility that does dummy computing for max 2 seconds to spin up cores
344 * This routine will check the number of cores configured and online
345 * (using sysconf), and the spins doing dummy compute operations for up to
346 * 2 seconds, or until all cores have come online. This can be used prior to
347 * hardware detection for platforms that take unused processors offline.
349 * This routine will not throw exceptions.
351 static void
352 spinUpCore() noexcept
354 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROCESSORS_ONLN)
355 float dummy = 0.1;
356 int countConfigured = sysconf(_SC_NPROCESSORS_CONF); // noexcept
357 auto start = std::chrono::steady_clock::now(); // noexcept
359 while (sysconf(_SC_NPROCESSORS_ONLN) < countConfigured &&
360 std::chrono::steady_clock::now() - start < std::chrono::seconds(2))
362 for (int i = 1; i < 10000; i++)
364 dummy /= i;
368 if (dummy < 0)
370 printf("This cannot happen, but prevents loop from being optimized away.");
372 #endif
375 /*! \brief Prepare the system before hardware topology detection
377 * This routine should perform any actions we want to put the system in a state
378 * where we want it to be before detecting the hardware topology. For most
379 * processors there is nothing to do, but some architectures (in particular ARM)
380 * have support for taking configured cores offline, which will make them disappear
381 * from the online processor count.
383 * This routine checks if there is a mismatch between the number of cores
384 * configured and online, and in that case we issue a small workload that
385 * attempts to wake sleeping cores before doing the actual detection.
387 * This type of mismatch can also occur for x86 or PowerPC on Linux, if SMT has only
388 * been disabled in the kernel (rather than bios). Since those cores will never
389 * come online automatically, we currently skip this test for x86 & PowerPC to
390 * avoid wasting 2 seconds. We also skip the test if there is no thread support.
392 * \note Cores will sleep relatively quickly again, so it's important to issue
393 * the real detection code directly after this routine.
395 static void
396 hardwareTopologyPrepareDetection()
398 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && \
399 (defined(THREAD_PTHREADS) || defined(THREAD_WINDOWS))
401 // Modify this conditional when/if x86 or PowerPC starts to sleep some cores
402 if (!isX86 && !isPowerPC)
404 int countConfigured = sysconf(_SC_NPROCESSORS_CONF);
405 std::vector<std::thread> workThreads(countConfigured);
407 for (auto &t : workThreads)
409 t = std::thread(spinUpCore);
412 for (auto &t : workThreads)
414 t.join();
417 #endif
420 /*! \brief Sanity check hardware topology and print some notes to log
422 * \param mdlog Logger.
423 * \param hardwareTopology Reference to hardwareTopology object.
425 static void
426 hardwareTopologyDoubleCheckDetection(const gmx::MDLogger gmx_unused &mdlog,
427 const gmx::HardwareTopology gmx_unused &hardwareTopology)
429 #if defined HAVE_SYSCONF && defined(_SC_NPROCESSORS_CONF)
430 if (hardwareTopology.supportLevel() < gmx::HardwareTopology::SupportLevel::LogicalProcessorCount)
432 return;
435 int countFromDetection = hardwareTopology.machine().logicalProcessorCount;
436 int countConfigured = sysconf(_SC_NPROCESSORS_CONF);
438 /* BIOS, kernel or user actions can take physical processors
439 * offline. We already cater for the some of the cases inside the hardwareToplogy
440 * by trying to spin up cores just before we detect, but there could be other
441 * cases where it is worthwhile to hint that there might be more resources available.
443 if (countConfigured >= 0 && countConfigured != countFromDetection)
445 GMX_LOG(mdlog.info).
446 appendTextFormatted("Note: %d CPUs configured, but only %d were detected to be online.\n", countConfigured, countFromDetection);
448 if (isX86 && countConfigured == 2*countFromDetection)
450 GMX_LOG(mdlog.info).
451 appendText(" X86 Hyperthreading is likely disabled; enable it for better performance.");
453 // For PowerPC (likely Power8) it is possible to set SMT to either 2,4, or 8-way hardware threads.
454 // We only warn if it is completely disabled since default performance drops with SMT8.
455 if (isPowerPC && countConfigured == 8*countFromDetection)
457 GMX_LOG(mdlog.info).
458 appendText(" PowerPC SMT is likely disabled; enable SMT2/SMT4 for better performance.");
461 #endif
465 gmx_hw_info_t *gmx_detect_hardware(const gmx::MDLogger &mdlog, const t_commrec *cr,
466 gmx_bool bDetectGPUs)
468 int ret;
470 /* make sure no one else is doing the same thing */
471 ret = tMPI_Thread_mutex_lock(&hw_info_lock);
472 if (ret != 0)
474 gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
477 /* only initialize the hwinfo structure if it is not already initalized */
478 if (n_hwinfo == 0)
480 snew(hwinfo_g, 1);
482 hwinfo_g->cpuInfo = new gmx::CpuInfo(gmx::CpuInfo::detect());
484 hardwareTopologyPrepareDetection();
485 hwinfo_g->hardwareTopology = new gmx::HardwareTopology(gmx::HardwareTopology::detect());
487 // If we detected the topology on this system, double-check that it makes sense
488 if (hwinfo_g->hardwareTopology->isThisSystem())
490 hardwareTopologyDoubleCheckDetection(mdlog, *(hwinfo_g->hardwareTopology));
493 // TODO: Get rid of this altogether.
494 hwinfo_g->nthreads_hw_avail = hwinfo_g->hardwareTopology->machine().logicalProcessorCount;
496 /* detect GPUs */
497 hwinfo_g->gpu_info.n_dev = 0;
498 hwinfo_g->gpu_info.n_dev_compatible = 0;
499 hwinfo_g->gpu_info.gpu_dev = nullptr;
501 /* Run the detection if the binary was compiled with GPU support
502 * and we requested detection.
504 hwinfo_g->gpu_info.bDetectGPUs =
505 (bGPUBinary && bDetectGPUs &&
506 getenv("GMX_DISABLE_GPU_DETECTION") == nullptr);
507 if (hwinfo_g->gpu_info.bDetectGPUs)
509 gmx_detect_gpus(mdlog, cr);
512 gmx_collect_hardware_mpi(*hwinfo_g->cpuInfo);
514 /* increase the reference counter */
515 n_hwinfo++;
517 ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
518 if (ret != 0)
520 gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
523 return hwinfo_g;
526 bool compatibleGpusFound(const gmx_gpu_info_t &gpu_info)
528 return gpu_info.n_dev_compatible > 0;
531 void gmx_hardware_info_free(gmx_hw_info_t *hwinfo)
533 int ret;
535 ret = tMPI_Thread_mutex_lock(&hw_info_lock);
536 if (ret != 0)
538 gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
541 /* decrease the reference counter */
542 n_hwinfo--;
545 if (hwinfo != hwinfo_g)
547 gmx_incons("hwinfo < hwinfo_g");
550 if (n_hwinfo < 0)
552 gmx_incons("n_hwinfo < 0");
555 if (n_hwinfo == 0)
557 delete hwinfo_g->cpuInfo;
558 delete hwinfo_g->hardwareTopology;
559 free_gpu_info(&hwinfo_g->gpu_info);
560 sfree(hwinfo_g);
563 ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
564 if (ret != 0)
566 gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));