Avoid GPU selection with GPU emulation
[alexxy/gromacs.git] / src / gromacs / gmxlib / gmx_detect_hardware.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015, 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 "gromacs/legacyheaders/gmx_detect_hardware.h"
38
39 #include "config.h"
40
41 #include <assert.h>
42 #include <errno.h>
43 #include <stdlib.h>
44 #include <string.h>
45
46 #include <string>
47 #include <vector>
48
49 #ifdef HAVE_UNISTD_H
50 /* For sysconf */
51 #include <unistd.h>
52 #endif
53 #ifdef GMX_NATIVE_WINDOWS
54 #include <windows.h>
55 #endif
56
57 #include "thread_mpi/threads.h"
58
59 #include "gromacs/gmxlib/gpu_utils/gpu_utils.h"
60 #include "gromacs/legacyheaders/copyrite.h"
61 #include "gromacs/legacyheaders/gmx_cpuid.h"
62 #include "gromacs/legacyheaders/md_logging.h"
63 #include "gromacs/legacyheaders/network.h"
64 #include "gromacs/legacyheaders/types/commrec.h"
65 #include "gromacs/legacyheaders/types/enums.h"
66 #include "gromacs/legacyheaders/types/hw_info.h"
67 #include "gromacs/utility/arrayref.h"
68 #include "gromacs/utility/basenetwork.h"
69 #include "gromacs/utility/cstringutil.h"
70 #include "gromacs/utility/exceptions.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/gmxassert.h"
73 #include "gromacs/utility/gmxomp.h"
74 #include "gromacs/utility/smalloc.h"
75 #include "gromacs/utility/stringutil.h"
76 #include "gromacs/utility/sysinfo.h"
77
78
79 #ifdef GMX_GPU
80 const gmx_bool bGPUBinary = TRUE;
81 #else
82 const gmx_bool bGPUBinary = FALSE;
83 #endif
84
85 static const char * invalid_gpuid_hint =
86     "A delimiter-free sequence of valid numeric IDs of available GPUs is expected.";
87
88 /* The globally shared hwinfo structure. */
89 static gmx_hw_info_t      *hwinfo_g;
90 /* A reference counter for the hwinfo structure */
91 static int                 n_hwinfo = 0;
92 /* A lock to protect the hwinfo structure */
93 static tMPI_Thread_mutex_t hw_info_lock = TMPI_THREAD_MUTEX_INITIALIZER;
94
95 #define HOSTNAMELEN 80
96
97 /* FW decl. */
98 static void set_gpu_ids(gmx_gpu_opt_t *gpu_opt, int nrank, int rank);
99 static int gmx_count_gpu_dev_unique(const gmx_gpu_info_t *gpu_info,
100                                     const gmx_gpu_opt_t  *gpu_opt);
101
102 static void sprint_gpus(char *sbuf, const gmx_gpu_info_t *gpu_info)
103 {
104     int      i, ndev;
105     char     stmp[STRLEN];
106
107     ndev = gpu_info->n_dev;
108
109     sbuf[0] = '\0';
110     for (i = 0; i < ndev; i++)
111     {
112         get_gpu_device_info_string(stmp, gpu_info, i);
113         strcat(sbuf, "    ");
114         strcat(sbuf, stmp);
115         if (i < ndev - 1)
116         {
117             strcat(sbuf, "\n");
118         }
119     }
120 }
121
122 static void print_gpu_detection_stats(FILE                 *fplog,
123                                       const gmx_gpu_info_t *gpu_info,
124                                       const t_commrec      *cr)
125 {
126     char onhost[HOSTNAMELEN+10], stmp[STRLEN];
127     int  ngpu;
128
129     if (!gpu_info->bDetectGPUs)
130     {
131         /* We skipped the detection, so don't print detection stats */
132         return;
133     }
134
135     ngpu = gpu_info->n_dev;
136
137 #if defined GMX_MPI && !defined GMX_THREAD_MPI
138     /* We only print the detection on one, of possibly multiple, nodes */
139     strncpy(onhost, " on host ", 10);
140     gmx_gethostname(onhost + 9, HOSTNAMELEN);
141 #else
142     /* We detect all relevant GPUs */
143     strncpy(onhost, "", 1);
144 #endif
145
146     if (ngpu > 0)
147     {
148         sprint_gpus(stmp, gpu_info);
149         md_print_warn(cr, fplog, "%d GPU%s detected%s:\n%s\n",
150                       ngpu, (ngpu > 1) ? "s" : "", onhost, stmp);
151     }
152     else
153     {
154         md_print_warn(cr, fplog, "No GPUs detected%s\n", onhost);
155     }
156 }
157
158 /*! \brief Helper function for reporting GPU usage information
159  * in the mdrun log file
160  *
161  * \param[in] gpu_info       Pointer to per-node GPU info struct
162  * \param[in] gpu_opt        Pointer to per-node GPU options struct
163  * \param[in] numPpRanks     Number of PP ranks per node
164  * \param[in] bPrintHostName Print the hostname in the usage information
165  * \return                   String to write to the log file
166  * \throws                   std::bad_alloc if out of memory */
167 static std::string
168 makeGpuUsageReport(const gmx_gpu_info_t *gpu_info,
169                    const gmx_gpu_opt_t  *gpu_opt,
170                    size_t                numPpRanks,
171                    bool                  bPrintHostName)
172 {
173     int  ngpu_use  = gpu_opt->n_dev_use;
174     int  ngpu_comp = gpu_info->n_dev_compatible;
175     char host[HOSTNAMELEN];
176
177     if (bPrintHostName)
178     {
179         gmx_gethostname(host, HOSTNAMELEN);
180     }
181
182     /* Issue a note if GPUs are available but not used */
183     if (ngpu_comp > 0 && ngpu_use < 1)
184     {
185         return gmx::formatString("%d compatible GPU%s detected in the system, but none will be used.\n"
186                                  "Consider trying GPU acceleration with the Verlet scheme!\n",
187                                  ngpu_comp, (ngpu_comp > 1) ? "s" : "");
188     }
189
190     std::string output;
191     if (!gpu_opt->bUserSet)
192     {
193         // gpu_opt->dev_compatible is only populated during auto-selection
194         std::string gpuIdsString =
195             formatAndJoin(gmx::constArrayRefFromArray(gpu_opt->dev_compatible,
196                                                       gpu_opt->n_dev_compatible),
197                           ",", gmx::StringFormatter("%d"));
198         bool bPluralGpus = gpu_opt->n_dev_compatible > 1;
199
200         if (bPrintHostName)
201         {
202             output += gmx::formatString("On host %s ", host);
203         }
204         output += gmx::formatString("%d compatible GPU%s %s present, with ID%s %s\n",
205                                     gpu_opt->n_dev_compatible,
206                                     bPluralGpus ? "s" : "",
207                                     bPluralGpus ? "are" : "is",
208                                     bPluralGpus ? "s" : "",
209                                     gpuIdsString.c_str());
210     }
211
212     {
213         std::vector<int>   gpuIdsInUse;
214         for (int i = 0; i < ngpu_use; i++)
215         {
216             gpuIdsInUse.push_back(get_cuda_gpu_device_id(gpu_info, gpu_opt, i));
217         }
218         std::string gpuIdsString =
219             formatAndJoin(gpuIdsInUse, ",", gmx::StringFormatter("%d"));
220         int         numGpusInUse = gmx_count_gpu_dev_unique(gpu_info, gpu_opt);
221         bool        bPluralGpus  = numGpusInUse > 1;
222
223         if (bPrintHostName)
224         {
225             output += gmx::formatString("On host %s ", host);
226         }
227         output += gmx::formatString("%d GPU%s %sselected for this run.\n"
228                                     "Mapping of GPU ID%s to the %d PP rank%s in this node: %s\n",
229                                     numGpusInUse, bPluralGpus ? "s" : "",
230                                     gpu_opt->bUserSet ? "user-" : "auto-",
231                                     bPluralGpus ? "s" : "",
232                                     numPpRanks,
233                                     (numPpRanks > 1) ? "s" : "",
234                                     gpuIdsString.c_str());
235     }
236
237     return output;
238 }
239
240 /* Give a suitable fatal error or warning if the build configuration
241    and runtime CPU do not match. */
242 static void
243 check_use_of_rdtscp_on_this_cpu(FILE                *fplog,
244                                 const t_commrec     *cr,
245                                 const gmx_hw_info_t *hwinfo)
246 {
247     gmx_bool bCpuHasRdtscp, bBinaryUsesRdtscp;
248 #ifdef HAVE_RDTSCP
249     bBinaryUsesRdtscp = TRUE;
250 #else
251     bBinaryUsesRdtscp = FALSE;
252 #endif
253
254     bCpuHasRdtscp = gmx_cpuid_feature(hwinfo->cpuid_info, GMX_CPUID_FEATURE_X86_RDTSCP);
255
256     if (!bCpuHasRdtscp && bBinaryUsesRdtscp)
257     {
258         gmx_fatal(FARGS, "The %s executable was compiled to use the rdtscp CPU instruction. "
259                   "However, this is not supported by the current hardware and continuing would lead to a crash. "
260                   "Please rebuild GROMACS with the GMX_USE_RDTSCP=OFF CMake option.",
261                   ShortProgram());
262     }
263
264     if (bCpuHasRdtscp && !bBinaryUsesRdtscp)
265     {
266         md_print_warn(cr, fplog, "The current CPU can measure timings more accurately than the code in\n"
267                       "%s was configured to use. This might affect your simulation\n"
268                       "speed as accurate timings are needed for load-balancing.\n"
269                       "Please consider rebuilding %s with the GMX_USE_RDTSCP=ON CMake option.\n",
270                       ShortProgram(), ShortProgram());
271     }
272 }
273
274 void gmx_check_hw_runconf_consistency(FILE                *fplog,
275                                       const gmx_hw_info_t *hwinfo,
276                                       const t_commrec     *cr,
277                                       const gmx_hw_opt_t  *hw_opt,
278                                       gmx_bool             bUseGPU)
279 {
280     int      npppn;
281     char     th_or_proc[STRLEN], th_or_proc_plural[STRLEN], pernode[STRLEN];
282     gmx_bool btMPI, bMPI, bMaxMpiThreadsSet, bNthreadsAuto, bEmulateGPU;
283
284     assert(hwinfo);
285     assert(cr);
286
287     /* Below we only do consistency checks for PP and GPUs,
288      * this is irrelevant for PME only nodes, so in that case we return
289      * here.
290      */
291     if (!(cr->duty & DUTY_PP))
292     {
293         return;
294     }
295
296 #if defined(GMX_THREAD_MPI)
297     bMPI          = FALSE;
298     btMPI         = TRUE;
299     bNthreadsAuto = (hw_opt->nthreads_tmpi < 1);
300 #elif defined(GMX_LIB_MPI)
301     bMPI          = TRUE;
302     btMPI         = FALSE;
303     bNthreadsAuto = FALSE;
304 #else
305     bMPI          = FALSE;
306     btMPI         = FALSE;
307     bNthreadsAuto = FALSE;
308 #endif
309
310     /* GPU emulation detection is done later, but we need here as well
311      * -- uncool, but there's no elegant workaround */
312     bEmulateGPU       = (getenv("GMX_EMULATE_GPU") != NULL);
313     bMaxMpiThreadsSet = (getenv("GMX_MAX_MPI_THREADS") != NULL);
314
315     if (hwinfo->gpu_info.n_dev_compatible > 0)
316     {
317         std::string gpuUseageReport;
318         try
319         {
320             gpuUseageReport = makeGpuUsageReport(&hwinfo->gpu_info,
321                                                  &hw_opt->gpu_opt,
322                                                  cr->nrank_pp_intranode,
323                                                  bMPI && cr->nnodes > 1);
324         }
325         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
326
327         /* NOTE: this print is only for and on one physical node */
328         md_print_info(cr, fplog, "%s\n", gpuUseageReport.c_str());
329     }
330
331     /* Need to ensure that we have enough GPUs:
332      * - need one GPU per PP node
333      * - no GPU oversubscription with tMPI
334      * */
335     /* number of PP processes per node */
336     npppn = cr->nrank_pp_intranode;
337
338     pernode[0]           = '\0';
339     th_or_proc_plural[0] = '\0';
340     if (btMPI)
341     {
342         sprintf(th_or_proc, "thread-MPI thread");
343         if (npppn > 1)
344         {
345             sprintf(th_or_proc_plural, "s");
346         }
347     }
348     else if (bMPI)
349     {
350         sprintf(th_or_proc, "MPI process");
351         if (npppn > 1)
352         {
353             sprintf(th_or_proc_plural, "es");
354         }
355         sprintf(pernode, " per node");
356     }
357     else
358     {
359         /* neither MPI nor tMPI */
360         sprintf(th_or_proc, "process");
361     }
362
363     if (bUseGPU && hwinfo->gpu_info.n_dev_compatible > 0 &&
364         !bEmulateGPU)
365     {
366         int  ngpu_comp, ngpu_use;
367         char gpu_comp_plural[2], gpu_use_plural[2];
368
369         ngpu_comp = hwinfo->gpu_info.n_dev_compatible;
370         ngpu_use  = hw_opt->gpu_opt.n_dev_use;
371
372         sprintf(gpu_comp_plural, "%s", (ngpu_comp > 1) ? "s" : "");
373         sprintf(gpu_use_plural,  "%s", (ngpu_use > 1) ? "s" : "");
374
375         /* number of tMPI threads auto-adjusted */
376         if (btMPI && bNthreadsAuto)
377         {
378             if (hw_opt->gpu_opt.bUserSet && npppn < ngpu_use)
379             {
380                 /* The user manually provided more GPUs than threads we
381                    could automatically start. */
382                 gmx_fatal(FARGS,
383                           "%d GPU%s provided, but only %d PP thread-MPI thread%s coud be started.\n"
384                           "%s requires one PP tread-MPI thread per GPU; use fewer GPUs%s.",
385                           ngpu_use, gpu_use_plural,
386                           npppn, th_or_proc_plural,
387                           ShortProgram(), bMaxMpiThreadsSet ? "\nor allow more threads to be used" : "");
388             }
389
390             if (!hw_opt->gpu_opt.bUserSet && npppn < ngpu_comp)
391             {
392                 /* There are more GPUs than tMPI threads; we have
393                    limited the number GPUs used. */
394                 md_print_warn(cr, fplog,
395                               "NOTE: %d GPU%s were detected, but only %d PP thread-MPI thread%s can be started.\n"
396                               "      %s can use one GPU per PP tread-MPI thread, so only %d GPU%s will be used.%s\n",
397                               ngpu_comp, gpu_comp_plural,
398                               npppn, th_or_proc_plural,
399                               ShortProgram(), npppn,
400                               npppn > 1 ? "s" : "",
401                               bMaxMpiThreadsSet ? "\n      Also, you can allow more threads to be used by increasing GMX_MAX_MPI_THREADS" : "");
402             }
403         }
404
405         if (hw_opt->gpu_opt.bUserSet)
406         {
407             if (ngpu_use != npppn)
408             {
409                 gmx_fatal(FARGS,
410                           "Incorrect launch configuration: mismatching number of PP %s%s and GPUs%s.\n"
411                           "%s was started with %d PP %s%s%s, but you provided %d GPU%s.",
412                           th_or_proc, btMPI ? "s" : "es", pernode,
413                           ShortProgram(), npppn, th_or_proc,
414                           th_or_proc_plural, pernode,
415                           ngpu_use, gpu_use_plural);
416             }
417         }
418         else
419         {
420             if (ngpu_comp > npppn)
421             {
422                 md_print_warn(cr, fplog,
423                               "NOTE: potentially sub-optimal launch configuration, %s started with less\n"
424                               "      PP %s%s%s than GPU%s available.\n"
425                               "      Each PP %s can use only one GPU, %d GPU%s%s will be used.\n",
426                               ShortProgram(), th_or_proc,
427                               th_or_proc_plural, pernode, gpu_comp_plural,
428                               th_or_proc, npppn, gpu_use_plural, pernode);
429             }
430
431             if (ngpu_use != npppn)
432             {
433                 /* Avoid duplicate error messages.
434                  * Unfortunately we can only do this at the physical node
435                  * level, since the hardware setup and MPI process count
436                  * might differ between physical nodes.
437                  */
438                 if (cr->rank_pp_intranode == 0)
439                 {
440                     gmx_fatal(FARGS,
441                               "Incorrect launch configuration: mismatching number of PP %s%s and GPUs%s.\n"
442                               "%s was started with %d PP %s%s%s, but only %d GPU%s were detected.",
443                               th_or_proc, btMPI ? "s" : "es", pernode,
444                               ShortProgram(), npppn, th_or_proc,
445                               th_or_proc_plural, pernode,
446                               ngpu_use, gpu_use_plural);
447                 }
448             }
449         }
450
451         {
452             int      same_count;
453
454             same_count = gmx_count_gpu_dev_shared(&hw_opt->gpu_opt);
455
456             if (same_count > 0)
457             {
458                 md_print_info(cr, fplog,
459                               "NOTE: You assigned %s to multiple %s%s.\n",
460                               same_count > 1 ? "GPUs" : "a GPU", th_or_proc, btMPI ? "s" : "es");
461             }
462         }
463     }
464
465 #ifdef GMX_MPI
466     if (PAR(cr))
467     {
468         /* Avoid other ranks to continue after
469            inconsistency */
470         MPI_Barrier(cr->mpi_comm_mygroup);
471     }
472 #endif
473
474 }
475
476 /* Return 0 if none of the GPU (per node) are shared among PP ranks.
477  *
478  * Sharing GPUs among multiple PP ranks is possible when the user passes
479  * GPU IDs. Here we check for sharing and return a non-zero value when
480  * this is detected. Note that the return value represents the number of
481  * PP rank pairs that share a device.
482  */
483 int gmx_count_gpu_dev_shared(const gmx_gpu_opt_t *gpu_opt)
484 {
485     int      same_count    = 0;
486     int      ngpu          = gpu_opt->n_dev_use;
487
488     if (gpu_opt->bUserSet)
489     {
490         int      i, j;
491
492         for (i = 0; i < ngpu - 1; i++)
493         {
494             for (j = i + 1; j < ngpu; j++)
495             {
496                 same_count      += (gpu_opt->dev_use[i] ==
497                                     gpu_opt->dev_use[j]);
498             }
499         }
500     }
501
502     return same_count;
503 }
504
505 /* Count and return the number of unique GPUs (per node) selected.
506  *
507  * As sharing GPUs among multiple PP ranks is possible when the user passes
508  * GPU IDs, the number of GPUs user (per node) can be different from the
509  * number of GPU IDs selected.
510  */
511 static int gmx_count_gpu_dev_unique(const gmx_gpu_info_t *gpu_info,
512                                     const gmx_gpu_opt_t  *gpu_opt)
513 {
514     int  i, uniq_count, ngpu;
515     int *uniq_ids;
516
517     assert(gpu_info);
518     assert(gpu_opt);
519
520     ngpu = gpu_info->n_dev;
521
522     uniq_count  = 0;
523
524     snew(uniq_ids, ngpu);
525
526     /* Each element in uniq_ids will be set to 0 or 1. The n-th element set
527      * to 1 indicates that the respective GPU was selected to be used. */
528     for (i = 0; i < gpu_opt->n_dev_use; i++)
529     {
530         uniq_ids[get_cuda_gpu_device_id(gpu_info, gpu_opt, i)] = 1;
531     }
532     /* Count the devices used. */
533     for (i = 0; i < ngpu; i++)
534     {
535         uniq_count += uniq_ids[i];
536     }
537
538     sfree(uniq_ids);
539
540     return uniq_count;
541 }
542
543 static int get_ncores(gmx_cpuid_t cpuid)
544 {
545     int        nprocessors, npackages, ncores_per_package, nhwthreads_per_core;
546     const int *package_id, *core_id, *hwthread_id, *locality_order;
547     int        rc;
548
549     rc = gmx_cpuid_topology(cpuid,
550                             &nprocessors, &npackages,
551                             &ncores_per_package, &nhwthreads_per_core,
552                             &package_id, &core_id,
553                             &hwthread_id, &locality_order);
554
555     if (rc == 0)
556     {
557         return npackages*ncores_per_package;
558     }
559     else
560     {
561         /* We don't have cpuid topology info, return 0 core count */
562         return 0;
563     }
564 }
565
566 /* Return the number of hardware threads supported by the current CPU.
567  * We assume that this is equal with the number of "processors"
568  * reported to be online by the OS at the time of the call. The
569  * definition of "processor" is according to an old POSIX standard.
570  *
571  * Note that the number of hardware threads is generally greater than
572  * the number of cores (e.g. x86 hyper-threading, Power). Managing the
573  * mapping of software threads to hardware threads is managed
574  * elsewhere. */
575 static int get_nthreads_hw_avail(FILE gmx_unused *fplog, const t_commrec gmx_unused *cr)
576 {
577     int ret = 0;
578
579 #if ((defined(WIN32) || defined( _WIN32 ) || defined(WIN64) || defined( _WIN64 )) && !(defined (__CYGWIN__) || defined (__CYGWIN32__)))
580     /* Windows */
581     SYSTEM_INFO sysinfo;
582     GetSystemInfo( &sysinfo );
583     ret = sysinfo.dwNumberOfProcessors;
584 #elif defined HAVE_SYSCONF
585     /* We are probably on Unix.
586      * Now check if we have the argument to use before executing the call
587      */
588 #if defined(_SC_NPROCESSORS_ONLN)
589     ret = sysconf(_SC_NPROCESSORS_ONLN);
590 #elif defined(_SC_NPROC_ONLN)
591     ret = sysconf(_SC_NPROC_ONLN);
592 #elif defined(_SC_NPROCESSORS_CONF)
593     ret = sysconf(_SC_NPROCESSORS_CONF);
594 #elif defined(_SC_NPROC_CONF)
595     ret = sysconf(_SC_NPROC_CONF);
596 #else
597 #warning "No valid sysconf argument value found. Executables will not be able to determine the number of logical cores: mdrun will use 1 thread by default!"
598 #endif /* End of check for sysconf argument values */
599
600 #else
601     /* Neither windows nor Unix. No fscking idea how many hardware threads we have! */
602     ret = -1;
603 #endif
604
605     if (debug)
606     {
607         fprintf(debug, "Detected %d hardware threads to use.\n", ret);
608     }
609
610 #ifdef GMX_OPENMP
611     if (ret != gmx_omp_get_num_procs())
612     {
613         md_print_warn(cr, fplog,
614                       "Number of logical cores detected (%d) does not match the number reported by OpenMP (%d).\n"
615                       "Consider setting the launch configuration manually!",
616                       ret, gmx_omp_get_num_procs());
617     }
618 #endif
619
620     return ret;
621 }
622
623 static void gmx_detect_gpus(FILE *fplog, const t_commrec *cr)
624 {
625 #ifdef GMX_LIB_MPI
626     int              rank_world;
627     MPI_Comm         physicalnode_comm;
628 #endif
629     int              rank_local;
630
631     /* Under certain circumstances MPI ranks on the same physical node
632      * can not simultaneously access the same GPU(s). Therefore we run
633      * the detection only on one MPI rank per node and broadcast the info.
634      * Note that with thread-MPI only a single thread runs this code.
635      *
636      * TODO: We should also do CPU hardware detection only once on each
637      * physical node and broadcast it, instead of do it on every MPI rank.
638      */
639 #ifdef GMX_LIB_MPI
640     /* A split of MPI_COMM_WORLD over physical nodes is only required here,
641      * so we create and destroy it locally.
642      */
643     MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
644     MPI_Comm_split(MPI_COMM_WORLD, gmx_physicalnode_id_hash(),
645                    rank_world, &physicalnode_comm);
646     MPI_Comm_rank(physicalnode_comm, &rank_local);
647 #else
648     /* Here there should be only one process, check this */
649     assert(cr->nnodes == 1 && cr->sim_nodeid == 0);
650
651     rank_local = 0;
652 #endif
653
654     if (rank_local == 0)
655     {
656         char detection_error[STRLEN] = "", sbuf[STRLEN];
657
658         if (detect_gpus(&hwinfo_g->gpu_info, detection_error) != 0)
659         {
660             if (detection_error[0] != '\0')
661             {
662                 sprintf(sbuf, ":\n      %s\n", detection_error);
663             }
664             else
665             {
666                 sprintf(sbuf, ".");
667             }
668             md_print_warn(cr, fplog,
669                           "NOTE: Error occurred during GPU detection%s"
670                           "      Can not use GPU acceleration, will fall back to CPU kernels.\n",
671                           sbuf);
672         }
673     }
674
675 #ifdef GMX_LIB_MPI
676     /* Broadcast the GPU info to the other ranks within this node */
677     MPI_Bcast(&hwinfo_g->gpu_info.n_dev, 1, MPI_INT, 0, physicalnode_comm);
678
679     if (hwinfo_g->gpu_info.n_dev > 0)
680     {
681         int dev_size;
682
683         dev_size = hwinfo_g->gpu_info.n_dev*sizeof_gpu_dev_info();
684
685         if (rank_local > 0)
686         {
687             hwinfo_g->gpu_info.gpu_dev =
688                 (struct gmx_device_info_t *)malloc(dev_size);
689         }
690         MPI_Bcast(hwinfo_g->gpu_info.gpu_dev, dev_size, MPI_BYTE,
691                   0, physicalnode_comm);
692         MPI_Bcast(&hwinfo_g->gpu_info.n_dev_compatible, 1, MPI_INT,
693                   0, physicalnode_comm);
694     }
695
696     MPI_Comm_free(&physicalnode_comm);
697 #endif
698 }
699
700 static void gmx_collect_hardware_mpi()
701 {
702 #ifdef GMX_LIB_MPI
703     int  rank_id;
704     int  nrank, rank, ncore, nhwthread, ngpu, i;
705     int  gpu_hash;
706     int *buf, *all;
707
708     rank_id   = gmx_physicalnode_id_hash();
709     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
710     MPI_Comm_size(MPI_COMM_WORLD, &nrank);
711     ncore     = hwinfo_g->ncore;
712     nhwthread = hwinfo_g->nthreads_hw_avail;
713     ngpu      = hwinfo_g->gpu_info.n_dev_compatible;
714     /* Create a unique hash of the GPU type(s) in this node */
715     gpu_hash  = 0;
716     /* Here it might be better to only loop over the compatible GPU, but we
717      * don't have that information available and it would also require
718      * removing the device ID from the device info string.
719      */
720     for (i = 0; i < hwinfo_g->gpu_info.n_dev; i++)
721     {
722         char stmp[STRLEN];
723
724         /* Since the device ID is incorporated in the hash, the order of
725          * the GPUs affects the hash. Also two identical GPUs won't give
726          * a gpu_hash of zero after XORing.
727          */
728         get_gpu_device_info_string(stmp, &hwinfo_g->gpu_info, i);
729         gpu_hash ^= gmx_string_fullhash_func(stmp, gmx_string_hash_init);
730     }
731
732     snew(buf, nrank);
733     snew(all, nrank);
734     buf[rank] = rank_id;
735
736     MPI_Allreduce(buf, all, nrank, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
737
738     gmx_bool bFound;
739     int      nnode0, ncore0, nhwthread0, ngpu0, r;
740
741     bFound     = FALSE;
742     ncore0     = 0;
743     nnode0     = 0;
744     nhwthread0 = 0;
745     ngpu0      = 0;
746     for (r = 0; r < nrank; r++)
747     {
748         if (all[r] == rank_id)
749         {
750             if (!bFound && r == rank)
751             {
752                 /* We are the first rank in this physical node */
753                 nnode0     = 1;
754                 ncore0     = ncore;
755                 nhwthread0 = nhwthread;
756                 ngpu0      = ngpu;
757             }
758             bFound = TRUE;
759         }
760     }
761
762     sfree(buf);
763     sfree(all);
764
765     int sum[4], maxmin[10];
766
767     {
768         int buf[4];
769
770         /* Sum values from only intra-rank 0 so we get the sum over all nodes */
771         buf[0] = nnode0;
772         buf[1] = ncore0;
773         buf[2] = nhwthread0;
774         buf[3] = ngpu0;
775
776         MPI_Allreduce(buf, sum, 4, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
777     }
778
779     {
780         int buf[10];
781
782         /* Store + and - values for all ranks,
783          * so we can get max+min with one MPI call.
784          */
785         buf[0] = ncore;
786         buf[1] = nhwthread;
787         buf[2] = ngpu;
788         buf[3] = gmx_cpuid_simd_suggest(hwinfo_g->cpuid_info);
789         buf[4] = gpu_hash;
790         buf[5] = -buf[0];
791         buf[6] = -buf[1];
792         buf[7] = -buf[2];
793         buf[8] = -buf[3];
794         buf[9] = -buf[4];
795
796         MPI_Allreduce(buf, maxmin, 10, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
797     }
798
799     hwinfo_g->nphysicalnode       = sum[0];
800     hwinfo_g->ncore_tot           = sum[1];
801     hwinfo_g->ncore_min           = -maxmin[5];
802     hwinfo_g->ncore_max           = maxmin[0];
803     hwinfo_g->nhwthread_tot       = sum[2];
804     hwinfo_g->nhwthread_min       = -maxmin[6];
805     hwinfo_g->nhwthread_max       = maxmin[1];
806     hwinfo_g->ngpu_compatible_tot = sum[3];
807     hwinfo_g->ngpu_compatible_min = -maxmin[7];
808     hwinfo_g->ngpu_compatible_max = maxmin[2];
809     hwinfo_g->simd_suggest_min    = static_cast<enum gmx_cpuid_simd>(-maxmin[8]);
810     hwinfo_g->simd_suggest_max    = static_cast<enum gmx_cpuid_simd>(maxmin[3]);
811     hwinfo_g->bIdenticalGPUs      = (maxmin[4] == -maxmin[9]);
812 #else
813     /* All ranks use the same pointer, protect it with a mutex */
814     tMPI_Thread_mutex_lock(&hw_info_lock);
815     hwinfo_g->nphysicalnode       = 1;
816     hwinfo_g->ncore_tot           = hwinfo_g->ncore;
817     hwinfo_g->ncore_min           = hwinfo_g->ncore;
818     hwinfo_g->ncore_max           = hwinfo_g->ncore;
819     hwinfo_g->nhwthread_tot       = hwinfo_g->nthreads_hw_avail;
820     hwinfo_g->nhwthread_min       = hwinfo_g->nthreads_hw_avail;
821     hwinfo_g->nhwthread_max       = hwinfo_g->nthreads_hw_avail;
822     hwinfo_g->ngpu_compatible_tot = hwinfo_g->gpu_info.n_dev_compatible;
823     hwinfo_g->ngpu_compatible_min = hwinfo_g->gpu_info.n_dev_compatible;
824     hwinfo_g->ngpu_compatible_max = hwinfo_g->gpu_info.n_dev_compatible;
825     hwinfo_g->simd_suggest_min    = gmx_cpuid_simd_suggest(hwinfo_g->cpuid_info);
826     hwinfo_g->simd_suggest_max    = gmx_cpuid_simd_suggest(hwinfo_g->cpuid_info);
827     hwinfo_g->bIdenticalGPUs      = TRUE;
828     tMPI_Thread_mutex_unlock(&hw_info_lock);
829 #endif
830 }
831
832 gmx_hw_info_t *gmx_detect_hardware(FILE *fplog, const t_commrec *cr,
833                                    gmx_bool bDetectGPUs)
834 {
835     int ret;
836
837     /* make sure no one else is doing the same thing */
838     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
839     if (ret != 0)
840     {
841         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
842     }
843
844     /* only initialize the hwinfo structure if it is not already initalized */
845     if (n_hwinfo == 0)
846     {
847         snew(hwinfo_g, 1);
848
849         /* detect CPUID info; no fuss, we don't detect system-wide
850          * -- sloppy, but that's it for now */
851         if (gmx_cpuid_init(&hwinfo_g->cpuid_info) != 0)
852         {
853             gmx_fatal_collective(FARGS, cr, NULL, "CPUID detection failed!");
854         }
855
856         /* get the number of cores, will be 0 when not detected */
857         hwinfo_g->ncore             = get_ncores(hwinfo_g->cpuid_info);
858
859         /* detect number of hardware threads */
860         hwinfo_g->nthreads_hw_avail = get_nthreads_hw_avail(fplog, cr);
861
862         /* detect GPUs */
863         hwinfo_g->gpu_info.n_dev            = 0;
864         hwinfo_g->gpu_info.n_dev_compatible = 0;
865         hwinfo_g->gpu_info.gpu_dev          = NULL;
866
867         /* Run the detection if the binary was compiled with GPU support
868          * and we requested detection.
869          */
870         hwinfo_g->gpu_info.bDetectGPUs =
871             (bGPUBinary && bDetectGPUs &&
872              getenv("GMX_DISABLE_GPU_DETECTION") == NULL);
873         if (hwinfo_g->gpu_info.bDetectGPUs)
874         {
875             gmx_detect_gpus(fplog, cr);
876         }
877     }
878     /* increase the reference counter */
879     n_hwinfo++;
880
881     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
882     if (ret != 0)
883     {
884         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
885     }
886
887     gmx_collect_hardware_mpi();
888
889     return hwinfo_g;
890 }
891
892 static std::string detected_hardware_string(const gmx_hw_info_t *hwinfo,
893                                             bool                 bFullCpuInfo)
894 {
895     std::string s;
896
897     s  = gmx::formatString("\n");
898     s += gmx::formatString("Running on %d node%s with total",
899                            hwinfo->nphysicalnode,
900                            hwinfo->nphysicalnode == 1 ? "" : "s");
901     if (hwinfo->ncore_tot > 0)
902     {
903         s += gmx::formatString(" %d cores,", hwinfo->ncore_tot);
904     }
905     s += gmx::formatString(" %d logical cores", hwinfo->nhwthread_tot);
906     if (hwinfo->gpu_info.bDetectGPUs)
907     {
908         s += gmx::formatString(", %d compatible GPU%s",
909                                hwinfo->ngpu_compatible_tot,
910                                hwinfo->ngpu_compatible_tot == 1 ? "" : "s");
911     }
912     else if (bGPUBinary)
913     {
914         s += gmx::formatString(" (GPU detection deactivated)");
915     }
916     s += gmx::formatString("\n");
917
918     if (hwinfo->nphysicalnode > 1)
919     {
920         /* Print per node hardware feature counts */
921         if (hwinfo->ncore_max > 0)
922         {
923             s += gmx::formatString("  Cores per node:           %2d", hwinfo->ncore_min);
924             if (hwinfo->ncore_max > hwinfo->ncore_min)
925             {
926                 s += gmx::formatString(" - %2d", hwinfo->ncore_max);
927             }
928             s += gmx::formatString("\n");
929         }
930         s += gmx::formatString("  Logical cores per node:   %2d", hwinfo->nhwthread_min);
931         if (hwinfo->nhwthread_max > hwinfo->nhwthread_min)
932         {
933             s += gmx::formatString(" - %2d", hwinfo->nhwthread_max);
934         }
935         s += gmx::formatString("\n");
936         if (bGPUBinary)
937         {
938             s += gmx::formatString("  Compatible GPUs per node: %2d",
939                                    hwinfo->ngpu_compatible_min);
940             if (hwinfo->ngpu_compatible_max > hwinfo->ngpu_compatible_min)
941             {
942                 s += gmx::formatString(" - %2d", hwinfo->ngpu_compatible_max);
943             }
944             s += gmx::formatString("\n");
945             if (hwinfo->ngpu_compatible_tot > 0)
946             {
947                 if (hwinfo->bIdenticalGPUs)
948                 {
949                     s += gmx::formatString("  All nodes have identical type(s) of GPUs\n");
950                 }
951                 else
952                 {
953                     /* This message will also appear with identical GPU types
954                      * when at least one node has no GPU.
955                      */
956                     s += gmx::formatString("  Different nodes have different type(s) and/or order of GPUs\n");
957                 }
958             }
959         }
960     }
961
962 #ifdef GMX_LIB_MPI
963     char host[HOSTNAMELEN];
964     int  rank;
965
966     gmx_gethostname(host, HOSTNAMELEN);
967     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
968
969     s += gmx::formatString("Hardware detected on host %s (the node of MPI rank %d):\n",
970                            host, rank);
971 #else
972     s += gmx::formatString("Hardware detected:\n");
973 #endif
974     s += gmx::formatString("  CPU info:\n");
975     if (bFullCpuInfo)
976     {
977         char buf[1024];
978
979         gmx_cpuid_formatstring(hwinfo->cpuid_info, buf, 1023);
980         buf[1023] = '\0';
981
982         s += gmx::formatString("%s", buf);
983     }
984     else
985     {
986         s += gmx::formatString("    Vendor: %s\n",
987                                gmx_cpuid_vendor_string[gmx_cpuid_vendor(hwinfo->cpuid_info)]);
988         s += gmx::formatString("    Brand:  %s\n",
989                                gmx_cpuid_brand(hwinfo->cpuid_info));
990     }
991     s += gmx::formatString("    SIMD instructions most likely to fit this hardware: %s",
992                            gmx_cpuid_simd_string[hwinfo->simd_suggest_min]);
993     if (hwinfo->simd_suggest_max > hwinfo->simd_suggest_min)
994     {
995         s += gmx::formatString(" - %s",
996                                gmx_cpuid_simd_string[hwinfo->simd_suggest_max]);
997     }
998     s += gmx::formatString("\n");
999     s += gmx::formatString("    SIMD instructions selected at GROMACS compile time: %s\n",
1000                            gmx_cpuid_simd_string[gmx_compiled_simd()]);
1001     if (bGPUBinary && (hwinfo->ngpu_compatible_tot > 0 ||
1002                        hwinfo->gpu_info.n_dev > 0))
1003     {
1004         s += gmx::formatString("  GPU info:\n");
1005         s += gmx::formatString("    Number of GPUs detected: %d\n",
1006                                hwinfo->gpu_info.n_dev);
1007         if (hwinfo->gpu_info.n_dev > 0)
1008         {
1009             char buf[STRLEN];
1010
1011             sprint_gpus(buf, &hwinfo->gpu_info);
1012             s += gmx::formatString("%s\n", buf);
1013         }
1014     }
1015
1016     return s;
1017 }
1018
1019 void gmx_print_detected_hardware(FILE *fplog, const t_commrec *cr,
1020                                  const gmx_hw_info_t *hwinfo)
1021 {
1022     if (fplog != NULL)
1023     {
1024         std::string detected;
1025
1026         detected = detected_hardware_string(hwinfo, TRUE);
1027
1028         fprintf(fplog, "%s\n", detected.c_str());
1029     }
1030
1031     if (MULTIMASTER(cr))
1032     {
1033         std::string detected;
1034
1035         detected = detected_hardware_string(hwinfo, FALSE);
1036
1037         fprintf(stderr, "%s\n", detected.c_str());
1038     }
1039
1040     /* Check the compiled SIMD instruction set against that of the node
1041      * with the lowest SIMD level support.
1042      */
1043     gmx_cpuid_simd_check(hwinfo->simd_suggest_min, fplog, MULTIMASTER(cr));
1044
1045     /* For RDTSCP we only check on our local node and skip the MPI reduction */
1046     check_use_of_rdtscp_on_this_cpu(fplog, cr, hwinfo);
1047 }
1048
1049 void gmx_parse_gpu_ids(gmx_gpu_opt_t *gpu_opt)
1050 {
1051     char *env;
1052
1053     if (gpu_opt->gpu_id != NULL && !bGPUBinary)
1054     {
1055         gmx_fatal(FARGS, "GPU ID string set, but %s was compiled without GPU support!", ShortProgram());
1056     }
1057
1058     env = getenv("GMX_GPU_ID");
1059     if (env != NULL && gpu_opt->gpu_id != NULL)
1060     {
1061         gmx_fatal(FARGS, "GMX_GPU_ID and -gpu_id can not be used at the same time");
1062     }
1063     if (env == NULL)
1064     {
1065         env = gpu_opt->gpu_id;
1066     }
1067
1068     /* parse GPU IDs if the user passed any */
1069     if (env != NULL)
1070     {
1071         /* Parse a "plain" GPU ID string which contains a sequence of
1072          * digits corresponding to GPU IDs; the order will indicate
1073          * the process/tMPI thread - GPU assignment. */
1074         parse_digits_from_plain_string(env,
1075                                        &gpu_opt->n_dev_use,
1076                                        &gpu_opt->dev_use);
1077
1078         if (gpu_opt->n_dev_use == 0)
1079         {
1080             gmx_fatal(FARGS, "Empty GPU ID string encountered.\n%s\n",
1081                       invalid_gpuid_hint);
1082         }
1083
1084         gpu_opt->bUserSet = TRUE;
1085     }
1086 }
1087
1088 void gmx_select_gpu_ids(FILE *fplog, const t_commrec *cr,
1089                         const gmx_gpu_info_t *gpu_info,
1090                         gmx_bool bForceUseGPU,
1091                         gmx_gpu_opt_t *gpu_opt)
1092 {
1093     int              i;
1094     char             sbuf[STRLEN], stmp[STRLEN];
1095
1096     /* Bail if binary is not compiled with GPU acceleration, but this is either
1097      * explicitly (-nb gpu) or implicitly (gpu ID passed) requested. */
1098     if (bForceUseGPU && !bGPUBinary)
1099     {
1100         gmx_fatal(FARGS, "GPU acceleration requested, but %s was compiled without GPU support!", ShortProgram());
1101     }
1102
1103     if (!(cr->duty & DUTY_PP))
1104     {
1105         /* Our rank is not doing PP, we don't use a GPU */
1106         return;
1107     }
1108
1109     if (gpu_opt->bUserSet)
1110     {
1111         /* Check the GPU IDs passed by the user.
1112          * (GPU IDs have been parsed by gmx_parse_gpu_ids before)
1113          */
1114         int *checkres;
1115         int  res;
1116
1117         snew(checkres, gpu_opt->n_dev_use);
1118
1119         res = check_selected_gpus(checkres, gpu_info, gpu_opt);
1120
1121         if (!res)
1122         {
1123             print_gpu_detection_stats(fplog, gpu_info, cr);
1124
1125             sprintf(sbuf, "Some of the requested GPUs do not exist, behave strangely, or are not compatible:\n");
1126             for (i = 0; i < gpu_opt->n_dev_use; i++)
1127             {
1128                 if (checkres[i] != egpuCompatible)
1129                 {
1130                     sprintf(stmp, "    GPU #%d: %s\n",
1131                             gpu_opt->dev_use[i],
1132                             gpu_detect_res_str[checkres[i]]);
1133                     strcat(sbuf, stmp);
1134                 }
1135             }
1136             gmx_fatal(FARGS, "%s", sbuf);
1137         }
1138
1139         sfree(checkres);
1140     }
1141     else if (getenv("GMX_EMULATE_GPU") == NULL)
1142     {
1143         pick_compatible_gpus(&hwinfo_g->gpu_info, gpu_opt);
1144         set_gpu_ids(gpu_opt, cr->nrank_pp_intranode, cr->rank_pp_intranode);
1145     }
1146
1147     /* If the user asked for a GPU, check whether we have a GPU */
1148     if (bForceUseGPU && gpu_info->n_dev_compatible == 0)
1149     {
1150         gmx_fatal(FARGS, "GPU acceleration requested, but no compatible GPUs were detected.");
1151     }
1152 }
1153
1154 /* Select the GPUs we will use. This is an operation local to each physical
1155  * node. If we have less MPI ranks than GPUs, we will waste some GPUs.
1156  * nrank and rank are the rank count and id for PP processes in our node.
1157  */
1158 static void set_gpu_ids(gmx_gpu_opt_t *gpu_opt, int nrank, int rank)
1159 {
1160     GMX_RELEASE_ASSERT(gpu_opt, "Invalid gpu_opt pointer passed");
1161     GMX_RELEASE_ASSERT(nrank >= 1,
1162                        gmx::formatString("Invalid limit (%d) for the number of GPUs (detected %d compatible GPUs)",
1163                                          rank, gpu_opt->n_dev_compatible).c_str());
1164
1165     if (gpu_opt->n_dev_compatible == 0)
1166     {
1167         char host[HOSTNAMELEN];
1168
1169         gmx_gethostname(host, HOSTNAMELEN);
1170         gmx_fatal(FARGS, "A GPU was requested on host %s, but no compatible GPUs were detected. All nodes with PP ranks need to have GPUs. If you intended to use GPU acceleration in a parallel run, you can either avoid using the nodes that don't have GPUs or place PME ranks on these nodes.", host);
1171     }
1172
1173     int nshare;
1174
1175     nshare = 1;
1176     if (nrank > gpu_opt->n_dev_compatible)
1177     {
1178         if (nrank % gpu_opt->n_dev_compatible == 0)
1179         {
1180             nshare = nrank/gpu_opt->n_dev_compatible;
1181         }
1182         else
1183         {
1184             if (rank == 0)
1185             {
1186                 gmx_fatal(FARGS, "The number of MPI ranks (%d) in a physical node is not a multiple of the number of GPUs (%d). Select a different number of MPI ranks or use the -gpu_id option to manually specify the GPU to be used.",
1187                           nrank, gpu_opt->n_dev_compatible);
1188             }
1189
1190 #ifdef GMX_MPI
1191             /* We use a global barrier to prevent ranks from continuing with
1192              * an invalid setup.
1193              */
1194             MPI_Barrier(MPI_COMM_WORLD);
1195 #endif
1196         }
1197     }
1198
1199     /* Here we will waste GPUs when nrank < gpu_opt->n_dev_compatible */
1200     gpu_opt->n_dev_use = std::min(gpu_opt->n_dev_compatible*nshare, nrank);
1201     snew(gpu_opt->dev_use, gpu_opt->n_dev_use);
1202     for (int i = 0; i != gpu_opt->n_dev_use; ++i)
1203     {
1204         /* TODO: improve this implementation: either sort GPUs or remove the weakest here */
1205         gpu_opt->dev_use[i] = gpu_opt->dev_compatible[i/nshare];
1206     }
1207 }
1208
1209 void gmx_hardware_info_free(gmx_hw_info_t *hwinfo)
1210 {
1211     int ret;
1212
1213     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
1214     if (ret != 0)
1215     {
1216         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
1217     }
1218
1219     /* decrease the reference counter */
1220     n_hwinfo--;
1221
1222
1223     if (hwinfo != hwinfo_g)
1224     {
1225         gmx_incons("hwinfo < hwinfo_g");
1226     }
1227
1228     if (n_hwinfo < 0)
1229     {
1230         gmx_incons("n_hwinfo < 0");
1231     }
1232
1233     if (n_hwinfo == 0)
1234     {
1235         gmx_cpuid_done(hwinfo_g->cpuid_info);
1236         free_gpu_info(&hwinfo_g->gpu_info);
1237         sfree(hwinfo_g);
1238     }
1239
1240     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
1241     if (ret != 0)
1242     {
1243         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
1244     }
1245 }