5a17e09256ee0c255d888168ec6f833ef4c6244f
[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, 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 <string>
48 #include <thread>
49 #include <vector>
50
51 #include "thread_mpi/threads.h"
52
53 #include "gromacs/gmxlib/network.h"
54 #include "gromacs/gpu_utils/gpu_utils.h"
55 #include "gromacs/hardware/cpuinfo.h"
56 #include "gromacs/hardware/gpu_hw_info.h"
57 #include "gromacs/hardware/hardwaretopology.h"
58 #include "gromacs/hardware/hw_info.h"
59 #include "gromacs/mdtypes/commrec.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/simd/support.h"
62 #include "gromacs/utility/arrayref.h"
63 #include "gromacs/utility/basedefinitions.h"
64 #include "gromacs/utility/basenetwork.h"
65 #include "gromacs/utility/baseversion.h"
66 #include "gromacs/utility/cstringutil.h"
67 #include "gromacs/utility/exceptions.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/logger.h"
71 #include "gromacs/utility/programcontext.h"
72 #include "gromacs/utility/smalloc.h"
73 #include "gromacs/utility/stringutil.h"
74 #include "gromacs/utility/sysinfo.h"
75
76 #ifdef HAVE_UNISTD_H
77 #    include <unistd.h>       // sysconf()
78 #endif
79
80 //! Convenience macro to help us avoid ifdefs each time we use sysconf
81 #if !defined(_SC_NPROCESSORS_ONLN) && defined(_SC_NPROC_ONLN)
82 #    define _SC_NPROCESSORS_ONLN _SC_NPROC_ONLN
83 #endif
84
85 //! Convenience macro to help us avoid ifdefs each time we use sysconf
86 #if !defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROC_CONF)
87 #    define _SC_NPROCESSORS_CONF _SC_NPROC_CONF
88 #endif
89
90 #if defined (__i386__) || defined (__x86_64__) || defined (_M_IX86) || defined (_M_X64)
91 //! Constant used to help minimize preprocessed code
92 static const bool isX86 = true;
93 #else
94 //! Constant used to help minimize preprocessed code
95 static const bool isX86 = false;
96 #endif
97
98 #if defined __powerpc__ || defined __ppc__ || defined __PPC__
99 static const bool isPowerPC = true;
100 #else
101 static const bool isPowerPC = false;
102 #endif
103
104 //! Constant used to help minimize preprocessed code
105 static const bool bGPUBinary     = GMX_GPU != GMX_GPU_NONE;
106
107 /* Note that some of the following arrays must match the "GPU support
108  * enumeration" in src/config.h.cmakein, so that GMX_GPU looks up an
109  * array entry. */
110
111 /* Both CUDA and OpenCL (on the supported/tested platforms) supports
112  * GPU device sharing.
113  */
114 static const bool gpuSharingSupport[] = { false, true, true };
115 static const bool bGpuSharingSupported = gpuSharingSupport[GMX_GPU];
116
117 /* Both CUDA and OpenCL (on the tested/supported platforms) supports everything.
118  */
119 static const bool multiGpuSupport[] = {
120     false, true, true
121 };
122 static const bool bMultiGpuPerNodeSupported = multiGpuSupport[GMX_GPU];
123
124 /* Names of the GPU detection/check results (see e_gpu_detect_res_t in hw_info.h). */
125 const char * const gpu_detect_res_str[egpuNR] =
126 {
127     "compatible", "inexistent", "incompatible", "insane"
128 };
129
130 static const char * invalid_gpuid_hint =
131     "A delimiter-free sequence of valid numeric IDs of available GPUs is expected.";
132
133 /* The globally shared hwinfo structure. */
134 static gmx_hw_info_t      *hwinfo_g;
135 /* A reference counter for the hwinfo structure */
136 static int                 n_hwinfo = 0;
137 /* A lock to protect the hwinfo structure */
138 static tMPI_Thread_mutex_t hw_info_lock = TMPI_THREAD_MUTEX_INITIALIZER;
139
140 #define HOSTNAMELEN 80
141
142 /* FW decl. */
143 static int gmx_count_gpu_dev_unique(const gmx_gpu_info_t *gpu_info,
144                                     const gmx_gpu_opt_t  *gpu_opt);
145
146 gmx_bool gmx_multiple_gpu_per_node_supported()
147 {
148     return bMultiGpuPerNodeSupported;
149 }
150
151 gmx_bool gmx_gpu_sharing_supported()
152 {
153     return bGpuSharingSupported;
154 }
155
156 std::string sprint_gpus(const gmx_gpu_info_t *gpu_info)
157 {
158     char                     stmp[STRLEN];
159     std::vector<std::string> gpuStrings;
160     for (int i = 0; i < gpu_info->n_dev; i++)
161     {
162         get_gpu_device_info_string(stmp, gpu_info, i);
163         gpuStrings.push_back(gmx::formatString("    %s", stmp));
164     }
165     return gmx::joinStrings(gpuStrings, "\n");
166 }
167
168 /*! \brief Helper function for reporting GPU usage information
169  * in the mdrun log file
170  *
171  * \param[in] gpu_info       Pointer to per-node GPU info struct
172  * \param[in] gpu_opt        Pointer to per-node GPU options struct
173  * \param[in] numPpRanks     Number of PP ranks per node
174  * \param[in] bPrintHostName Print the hostname in the usage information
175  * \return                   String to write to the log file
176  * \throws                   std::bad_alloc if out of memory */
177 static std::string
178 makeGpuUsageReport(const gmx_gpu_info_t *gpu_info,
179                    const gmx_gpu_opt_t  *gpu_opt,
180                    size_t                numPpRanks,
181                    bool                  bPrintHostName)
182 {
183     int  ngpu_use  = gpu_opt->n_dev_use;
184     int  ngpu_comp = gpu_info->n_dev_compatible;
185     char host[HOSTNAMELEN];
186
187     if (bPrintHostName)
188     {
189         gmx_gethostname(host, HOSTNAMELEN);
190     }
191
192     /* Issue a note if GPUs are available but not used */
193     if (ngpu_comp > 0 && ngpu_use < 1)
194     {
195         return gmx::formatString("%d compatible GPU%s detected in the system, but none will be used.\n"
196                                  "Consider trying GPU acceleration with the Verlet scheme!\n",
197                                  ngpu_comp, (ngpu_comp > 1) ? "s" : "");
198     }
199
200     std::string output;
201     if (!gpu_opt->bUserSet)
202     {
203         // gpu_opt->dev_compatible is only populated during auto-selection
204         std::string gpuIdsString =
205             formatAndJoin(gmx::constArrayRefFromArray(gpu_opt->dev_compatible,
206                                                       gpu_opt->n_dev_compatible),
207                           ",", gmx::StringFormatter("%d"));
208         bool bPluralGpus = gpu_opt->n_dev_compatible > 1;
209
210         if (bPrintHostName)
211         {
212             output += gmx::formatString("On host %s ", host);
213         }
214         output += gmx::formatString("%d compatible GPU%s %s present, with ID%s %s\n",
215                                     gpu_opt->n_dev_compatible,
216                                     bPluralGpus ? "s" : "",
217                                     bPluralGpus ? "are" : "is",
218                                     bPluralGpus ? "s" : "",
219                                     gpuIdsString.c_str());
220     }
221
222     {
223         std::vector<int> gpuIdsInUse;
224         for (int i = 0; i < ngpu_use; i++)
225         {
226             gpuIdsInUse.push_back(get_gpu_device_id(gpu_info, gpu_opt, i));
227         }
228         std::string gpuIdsString =
229             formatAndJoin(gpuIdsInUse, ",", gmx::StringFormatter("%d"));
230         int         numGpusInUse = gmx_count_gpu_dev_unique(gpu_info, gpu_opt);
231         bool        bPluralGpus  = numGpusInUse > 1;
232
233         if (bPrintHostName)
234         {
235             output += gmx::formatString("On host %s ", host);
236         }
237         output += gmx::formatString("%d GPU%s %sselected for this run.\n"
238                                     "Mapping of GPU ID%s to the %d PP rank%s in this node: %s\n",
239                                     numGpusInUse, bPluralGpus ? "s" : "",
240                                     gpu_opt->bUserSet ? "user-" : "auto-",
241                                     bPluralGpus ? "s" : "",
242                                     numPpRanks,
243                                     (numPpRanks > 1) ? "s" : "",
244                                     gpuIdsString.c_str());
245     }
246
247     return output;
248 }
249
250 /* Give a suitable fatal error or warning if the build configuration
251    and runtime CPU do not match. */
252 static void
253 check_use_of_rdtscp_on_this_cpu(const gmx::MDLogger   &mdlog,
254                                 const gmx::CpuInfo    &cpuInfo)
255 {
256 #ifdef HAVE_RDTSCP
257     bool binaryUsesRdtscp = TRUE;
258 #else
259     bool binaryUsesRdtscp = FALSE;
260 #endif
261
262     const char *programName = gmx::getProgramContext().displayName();
263
264     if (cpuInfo.supportLevel() < gmx::CpuInfo::SupportLevel::Features)
265     {
266         if (binaryUsesRdtscp)
267         {
268             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
269                     "The %s executable was compiled to use the rdtscp CPU instruction. "
270                     "We cannot detect the features of your current CPU, but will proceed anyway. "
271                     "If you get a crash, rebuild GROMACS with the GMX_USE_RDTSCP=OFF CMake option.",
272                     programName);
273         }
274     }
275     else
276     {
277         bool cpuHasRdtscp = cpuInfo.feature(gmx::CpuInfo::Feature::X86_Rdtscp);
278
279         if (!cpuHasRdtscp && binaryUsesRdtscp)
280         {
281             gmx_fatal(FARGS, "The %s executable was compiled to use the rdtscp CPU instruction. "
282                       "However, this is not supported by the current hardware and continuing would lead to a crash. "
283                       "Please rebuild GROMACS with the GMX_USE_RDTSCP=OFF CMake option.",
284                       programName);
285         }
286
287         if (cpuHasRdtscp && !binaryUsesRdtscp)
288         {
289             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
290                     "The current CPU can measure timings more accurately than the code in\n"
291                     "%s was configured to use. This might affect your simulation\n"
292                     "speed as accurate timings are needed for load-balancing.\n"
293                     "Please consider rebuilding %s with the GMX_USE_RDTSCP=ON CMake option.",
294                     programName, programName);
295         }
296     }
297 }
298
299 void gmx_check_hw_runconf_consistency(const gmx::MDLogger &mdlog,
300                                       const gmx_hw_info_t *hwinfo,
301                                       const t_commrec     *cr,
302                                       const gmx_hw_opt_t  *hw_opt,
303                                       gmx_bool             bUseGPU)
304 {
305     int      npppn;
306     char     th_or_proc[STRLEN], th_or_proc_plural[STRLEN], pernode[STRLEN];
307     gmx_bool btMPI, bMPI, bNthreadsAuto, bEmulateGPU;
308
309     GMX_RELEASE_ASSERT(hwinfo, "hwinfo must be a non-NULL pointer");
310     GMX_RELEASE_ASSERT(cr, "cr must be a non-NULL pointer");
311
312     /* Below we only do consistency checks for PP and GPUs,
313      * this is irrelevant for PME only nodes, so in that case we return
314      * here.
315      */
316     if (!(cr->duty & DUTY_PP))
317     {
318         return;
319     }
320
321 #if GMX_THREAD_MPI
322     bMPI          = FALSE;
323     btMPI         = TRUE;
324     bNthreadsAuto = (hw_opt->nthreads_tmpi < 1);
325 #elif GMX_LIB_MPI
326     bMPI          = TRUE;
327     btMPI         = FALSE;
328     bNthreadsAuto = FALSE;
329 #else
330     bMPI          = FALSE;
331     btMPI         = FALSE;
332     bNthreadsAuto = FALSE;
333 #endif
334
335     /* GPU emulation detection is done later, but we need here as well
336      * -- uncool, but there's no elegant workaround */
337     bEmulateGPU       = (getenv("GMX_EMULATE_GPU") != NULL);
338
339     if (hwinfo->gpu_info.n_dev_compatible > 0)
340     {
341         std::string gpuUsageReport;
342         try
343         {
344             gpuUsageReport = makeGpuUsageReport(&hwinfo->gpu_info,
345                                                 &hw_opt->gpu_opt,
346                                                 cr->nrank_pp_intranode,
347                                                 bMPI && cr->nnodes > 1);
348         }
349         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
350
351         /* NOTE: this print is only for and on one physical node */
352         GMX_LOG(mdlog.warning).appendText(gpuUsageReport);
353     }
354
355     /* Need to ensure that we have enough GPUs:
356      * - need one GPU per PP node
357      * - no GPU oversubscription with tMPI
358      * */
359     /* number of PP processes per node */
360     npppn = cr->nrank_pp_intranode;
361
362     pernode[0]           = '\0';
363     th_or_proc_plural[0] = '\0';
364     if (btMPI)
365     {
366         sprintf(th_or_proc, "thread-MPI thread");
367         if (npppn > 1)
368         {
369             sprintf(th_or_proc_plural, "s");
370         }
371     }
372     else if (bMPI)
373     {
374         sprintf(th_or_proc, "MPI process");
375         if (npppn > 1)
376         {
377             sprintf(th_or_proc_plural, "es");
378         }
379         sprintf(pernode, " per node");
380     }
381     else
382     {
383         /* neither MPI nor tMPI */
384         sprintf(th_or_proc, "process");
385     }
386
387     if (bUseGPU && hwinfo->gpu_info.n_dev_compatible > 0 &&
388         !bEmulateGPU)
389     {
390         int  ngpu_comp, ngpu_use;
391         char gpu_comp_plural[2], gpu_use_plural[2];
392
393         ngpu_comp = hwinfo->gpu_info.n_dev_compatible;
394         ngpu_use  = hw_opt->gpu_opt.n_dev_use;
395
396         sprintf(gpu_comp_plural, "%s", (ngpu_comp > 1) ? "s" : "");
397         sprintf(gpu_use_plural,  "%s", (ngpu_use > 1) ? "s" : "");
398
399         const char *programName = gmx::getProgramContext().displayName();
400
401         /* number of tMPI threads auto-adjusted */
402         if (btMPI && bNthreadsAuto)
403         {
404             if (hw_opt->gpu_opt.bUserSet && npppn < ngpu_use)
405             {
406                 /* The user manually provided more GPUs than threads we
407                    could automatically start. */
408                 gmx_fatal(FARGS,
409                           "%d GPU%s provided, but only %d PP thread-MPI thread%s coud be started.\n"
410                           "%s requires one PP tread-MPI thread per GPU; use fewer GPUs.",
411                           ngpu_use, gpu_use_plural,
412                           npppn, th_or_proc_plural,
413                           programName);
414             }
415
416             if (!hw_opt->gpu_opt.bUserSet && npppn < ngpu_comp)
417             {
418                 /* There are more GPUs than tMPI threads; we have
419                    limited the number GPUs used. */
420                 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
421                         "NOTE: %d GPU%s were detected, but only %d PP thread-MPI thread%s can be started.\n"
422                         "      %s can use one GPU per PP tread-MPI thread, so only %d GPU%s will be used.",
423                         ngpu_comp, gpu_comp_plural,
424                         npppn, th_or_proc_plural,
425                         programName, npppn,
426                         npppn > 1 ? "s" : "");
427             }
428         }
429
430         if (hw_opt->gpu_opt.bUserSet)
431         {
432             if (ngpu_use != npppn)
433             {
434                 gmx_fatal(FARGS,
435                           "Incorrect launch configuration: mismatching number of PP %s%s and GPUs%s.\n"
436                           "%s was started with %d PP %s%s%s, but you provided %d GPU%s.",
437                           th_or_proc, btMPI ? "s" : "es", pernode,
438                           programName, npppn, th_or_proc,
439                           th_or_proc_plural, pernode,
440                           ngpu_use, gpu_use_plural);
441             }
442         }
443         else
444         {
445             /* TODO Should we have a gpu_opt->n_dev_supported field? */
446             if (ngpu_comp > npppn && gmx_multiple_gpu_per_node_supported())
447             {
448                 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
449                         "NOTE: potentially sub-optimal launch configuration, %s started with less\n"
450                         "      PP %s%s%s than GPU%s available.\n"
451                         "      Each PP %s can use only one GPU, %d GPU%s%s will be used.",
452                         programName, th_or_proc,
453                         th_or_proc_plural, pernode, gpu_comp_plural,
454                         th_or_proc, npppn, gpu_use_plural, pernode);
455             }
456
457             if (ngpu_use != npppn)
458             {
459                 /* Avoid duplicate error messages.
460                  * Unfortunately we can only do this at the physical node
461                  * level, since the hardware setup and MPI process count
462                  * might differ between physical nodes.
463                  */
464                 if (cr->rank_pp_intranode == 0)
465                 {
466                     std::string reasonForLimit;
467                     if (ngpu_comp > 1 &&
468                         ngpu_use == 1 &&
469                         !gmx_multiple_gpu_per_node_supported())
470                     {
471                         reasonForLimit  = "can be used by ";
472                         reasonForLimit += getGpuImplementationString();
473                         reasonForLimit += " in GROMACS";
474                     }
475                     else
476                     {
477                         reasonForLimit = "was detected";
478                     }
479                     gmx_fatal(FARGS,
480                               "Incorrect launch configuration: mismatching number of PP %s%s and GPUs%s.\n"
481                               "%s was started with %d PP %s%s%s, but only %d GPU%s %s.",
482                               th_or_proc, btMPI ? "s" : "es", pernode,
483                               programName, npppn, th_or_proc,
484                               th_or_proc_plural, pernode,
485                               ngpu_use, gpu_use_plural, reasonForLimit.c_str());
486                 }
487             }
488         }
489
490         {
491             int      same_count;
492
493             same_count = gmx_count_gpu_dev_shared(&hw_opt->gpu_opt);
494
495             if (same_count > 0)
496             {
497                 GMX_LOG(mdlog.warning).appendTextFormatted(
498                         "NOTE: You assigned %s to multiple %s%s.",
499                         same_count > 1 ? "GPUs" : "a GPU", th_or_proc, btMPI ? "s" : "es");
500             }
501         }
502     }
503
504 #if GMX_MPI
505     if (PAR(cr))
506     {
507         /* Avoid other ranks to continue after
508            inconsistency */
509         MPI_Barrier(cr->mpi_comm_mygroup);
510     }
511 #endif
512
513 }
514
515 /* Return 0 if none of the GPU (per node) are shared among PP ranks.
516  *
517  * Sharing GPUs among multiple PP ranks is possible when the user passes
518  * GPU IDs. Here we check for sharing and return a non-zero value when
519  * this is detected. Note that the return value represents the number of
520  * PP rank pairs that share a device.
521  */
522 int gmx_count_gpu_dev_shared(const gmx_gpu_opt_t *gpu_opt)
523 {
524     int      same_count    = 0;
525     int      ngpu          = gpu_opt->n_dev_use;
526
527     if (gpu_opt->bUserSet)
528     {
529         int      i, j;
530
531         for (i = 0; i < ngpu - 1; i++)
532         {
533             for (j = i + 1; j < ngpu; j++)
534             {
535                 same_count      += (gpu_opt->dev_use[i] ==
536                                     gpu_opt->dev_use[j]);
537             }
538         }
539     }
540
541     return same_count;
542 }
543
544 /* Count and return the number of unique GPUs (per node) selected.
545  *
546  * As sharing GPUs among multiple PP ranks is possible when the user passes
547  * GPU IDs, the number of GPUs user (per node) can be different from the
548  * number of GPU IDs selected.
549  */
550 static int gmx_count_gpu_dev_unique(const gmx_gpu_info_t *gpu_info,
551                                     const gmx_gpu_opt_t  *gpu_opt)
552 {
553     int  i, uniq_count, ngpu;
554     int *uniq_ids;
555
556     GMX_RELEASE_ASSERT(gpu_info, "gpu_info must be a non-NULL pointer");
557     GMX_RELEASE_ASSERT(gpu_opt, "gpu_opt must be a non-NULL pointer");
558
559     ngpu = gpu_info->n_dev;
560
561     uniq_count  = 0;
562
563     snew(uniq_ids, ngpu);
564
565     /* Each element in uniq_ids will be set to 0 or 1. The n-th element set
566      * to 1 indicates that the respective GPU was selected to be used. */
567     for (i = 0; i < gpu_opt->n_dev_use; i++)
568     {
569         int device_id;
570
571         device_id           = gmx_gpu_sharing_supported() ? get_gpu_device_id(gpu_info, gpu_opt, i) : i;
572         uniq_ids[device_id] = 1;
573     }
574     /* Count the devices used. */
575     for (i = 0; i < ngpu; i++)
576     {
577         uniq_count += uniq_ids[i];
578     }
579
580     sfree(uniq_ids);
581
582     return uniq_count;
583 }
584
585 static void gmx_detect_gpus(const gmx::MDLogger &mdlog, const t_commrec *cr)
586 {
587 #if GMX_LIB_MPI
588     int              rank_world;
589     MPI_Comm         physicalnode_comm;
590 #endif
591     int              rank_local;
592
593     /* Under certain circumstances MPI ranks on the same physical node
594      * can not simultaneously access the same GPU(s). Therefore we run
595      * the detection only on one MPI rank per node and broadcast the info.
596      * Note that with thread-MPI only a single thread runs this code.
597      *
598      * NOTE: We can't broadcast gpu_info with OpenCL as the device and platform
599      * ID stored in the structure are unique for each rank (even if a device
600      * is shared by multiple ranks).
601      *
602      * TODO: We should also do CPU hardware detection only once on each
603      * physical node and broadcast it, instead of do it on every MPI rank.
604      */
605 #if GMX_LIB_MPI
606     /* A split of MPI_COMM_WORLD over physical nodes is only required here,
607      * so we create and destroy it locally.
608      */
609     MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
610     MPI_Comm_split(MPI_COMM_WORLD, gmx_physicalnode_id_hash(),
611                    rank_world, &physicalnode_comm);
612     MPI_Comm_rank(physicalnode_comm, &rank_local);
613     GMX_UNUSED_VALUE(cr);
614 #else
615     /* Here there should be only one process, check this */
616     GMX_RELEASE_ASSERT(cr->nnodes == 1 && cr->sim_nodeid == 0, "Only a single (master) process should execute here");
617
618     rank_local = 0;
619 #endif
620
621     /*  With CUDA detect only on one rank per host, with OpenCL need do
622      *  the detection on all PP ranks */
623     bool isOpenclPpRank = ((GMX_GPU == GMX_GPU_OPENCL) && (cr->duty & DUTY_PP));
624
625     if (rank_local == 0 || isOpenclPpRank)
626     {
627         char detection_error[STRLEN] = "", sbuf[STRLEN];
628
629         if (detect_gpus(&hwinfo_g->gpu_info, detection_error) != 0)
630         {
631             if (detection_error[0] != '\0')
632             {
633                 sprintf(sbuf, ":\n      %s\n", detection_error);
634             }
635             else
636             {
637                 sprintf(sbuf, ".");
638             }
639             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
640                     "NOTE: Error occurred during GPU detection%s"
641                     "      Can not use GPU acceleration, will fall back to CPU kernels.",
642                     sbuf);
643         }
644     }
645
646 #if GMX_LIB_MPI
647     if (!isOpenclPpRank)
648     {
649         /* Broadcast the GPU info to the other ranks within this node */
650         MPI_Bcast(&hwinfo_g->gpu_info.n_dev, 1, MPI_INT, 0, physicalnode_comm);
651
652         if (hwinfo_g->gpu_info.n_dev > 0)
653         {
654             int dev_size;
655
656             dev_size = hwinfo_g->gpu_info.n_dev*sizeof_gpu_dev_info();
657
658             if (rank_local > 0)
659             {
660                 hwinfo_g->gpu_info.gpu_dev =
661                     (struct gmx_device_info_t *)malloc(dev_size);
662             }
663             MPI_Bcast(hwinfo_g->gpu_info.gpu_dev, dev_size, MPI_BYTE,
664                       0, physicalnode_comm);
665             MPI_Bcast(&hwinfo_g->gpu_info.n_dev_compatible, 1, MPI_INT,
666                       0, physicalnode_comm);
667         }
668     }
669
670     MPI_Comm_free(&physicalnode_comm);
671 #endif
672 }
673
674 static void gmx_collect_hardware_mpi(const gmx::CpuInfo &cpuInfo)
675 {
676     const int ncore = hwinfo_g->hardwareTopology->numberOfCores();
677 #if GMX_LIB_MPI
678     int       rank_id;
679     int       nrank, rank, nhwthread, ngpu, i;
680     int       gpu_hash;
681     int      *buf, *all;
682
683     rank_id   = gmx_physicalnode_id_hash();
684     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
685     MPI_Comm_size(MPI_COMM_WORLD, &nrank);
686     nhwthread = hwinfo_g->nthreads_hw_avail;
687     ngpu      = hwinfo_g->gpu_info.n_dev_compatible;
688     /* Create a unique hash of the GPU type(s) in this node */
689     gpu_hash  = 0;
690     /* Here it might be better to only loop over the compatible GPU, but we
691      * don't have that information available and it would also require
692      * removing the device ID from the device info string.
693      */
694     for (i = 0; i < hwinfo_g->gpu_info.n_dev; i++)
695     {
696         char stmp[STRLEN];
697
698         /* Since the device ID is incorporated in the hash, the order of
699          * the GPUs affects the hash. Also two identical GPUs won't give
700          * a gpu_hash of zero after XORing.
701          */
702         get_gpu_device_info_string(stmp, &hwinfo_g->gpu_info, i);
703         gpu_hash ^= gmx_string_fullhash_func(stmp, gmx_string_hash_init);
704     }
705
706     snew(buf, nrank);
707     snew(all, nrank);
708     buf[rank] = rank_id;
709
710     MPI_Allreduce(buf, all, nrank, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
711
712     gmx_bool bFound;
713     int      nnode0, ncore0, nhwthread0, ngpu0, r;
714
715     bFound     = FALSE;
716     ncore0     = 0;
717     nnode0     = 0;
718     nhwthread0 = 0;
719     ngpu0      = 0;
720     for (r = 0; r < nrank; r++)
721     {
722         if (all[r] == rank_id)
723         {
724             if (!bFound && r == rank)
725             {
726                 /* We are the first rank in this physical node */
727                 nnode0     = 1;
728                 ncore0     = ncore;
729                 nhwthread0 = nhwthread;
730                 ngpu0      = ngpu;
731             }
732             bFound = TRUE;
733         }
734     }
735
736     sfree(buf);
737     sfree(all);
738
739     int sum[4], maxmin[10];
740
741     {
742         int buf[4];
743
744         /* Sum values from only intra-rank 0 so we get the sum over all nodes */
745         buf[0] = nnode0;
746         buf[1] = ncore0;
747         buf[2] = nhwthread0;
748         buf[3] = ngpu0;
749
750         MPI_Allreduce(buf, sum, 4, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
751     }
752
753     {
754         int buf[10];
755
756         /* Store + and - values for all ranks,
757          * so we can get max+min with one MPI call.
758          */
759         buf[0] = ncore;
760         buf[1] = nhwthread;
761         buf[2] = ngpu;
762         buf[3] = static_cast<int>(gmx::simdSuggested(cpuInfo));
763         buf[4] = gpu_hash;
764         buf[5] = -buf[0];
765         buf[6] = -buf[1];
766         buf[7] = -buf[2];
767         buf[8] = -buf[3];
768         buf[9] = -buf[4];
769
770         MPI_Allreduce(buf, maxmin, 10, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
771     }
772
773     hwinfo_g->nphysicalnode       = sum[0];
774     hwinfo_g->ncore_tot           = sum[1];
775     hwinfo_g->ncore_min           = -maxmin[5];
776     hwinfo_g->ncore_max           = maxmin[0];
777     hwinfo_g->nhwthread_tot       = sum[2];
778     hwinfo_g->nhwthread_min       = -maxmin[6];
779     hwinfo_g->nhwthread_max       = maxmin[1];
780     hwinfo_g->ngpu_compatible_tot = sum[3];
781     hwinfo_g->ngpu_compatible_min = -maxmin[7];
782     hwinfo_g->ngpu_compatible_max = maxmin[2];
783     hwinfo_g->simd_suggest_min    = -maxmin[8];
784     hwinfo_g->simd_suggest_max    = maxmin[3];
785     hwinfo_g->bIdenticalGPUs      = (maxmin[4] == -maxmin[9]);
786 #else
787     /* All ranks use the same pointer, protected by a mutex in the caller */
788     hwinfo_g->nphysicalnode       = 1;
789     hwinfo_g->ncore_tot           = ncore;
790     hwinfo_g->ncore_min           = ncore;
791     hwinfo_g->ncore_max           = ncore;
792     hwinfo_g->nhwthread_tot       = hwinfo_g->nthreads_hw_avail;
793     hwinfo_g->nhwthread_min       = hwinfo_g->nthreads_hw_avail;
794     hwinfo_g->nhwthread_max       = hwinfo_g->nthreads_hw_avail;
795     hwinfo_g->ngpu_compatible_tot = hwinfo_g->gpu_info.n_dev_compatible;
796     hwinfo_g->ngpu_compatible_min = hwinfo_g->gpu_info.n_dev_compatible;
797     hwinfo_g->ngpu_compatible_max = hwinfo_g->gpu_info.n_dev_compatible;
798     hwinfo_g->simd_suggest_min    = static_cast<int>(simdSuggested(cpuInfo));
799     hwinfo_g->simd_suggest_max    = static_cast<int>(simdSuggested(cpuInfo));
800     hwinfo_g->bIdenticalGPUs      = TRUE;
801 #endif
802 }
803
804 /*! \brief Utility that does dummy computing for max 2 seconds to spin up cores
805  *
806  *  This routine will check the number of cores configured and online
807  *  (using sysconf), and the spins doing dummy compute operations for up to
808  *  2 seconds, or until all cores have come online. This can be used prior to
809  *  hardware detection for platforms that take unused processors offline.
810  *
811  *  This routine will not throw exceptions.
812  */
813 static void
814 spinUpCore() noexcept
815 {
816 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && defined(_SC_NPROCESSORS_ONLN)
817     float dummy           = 0.1;
818     int   countConfigured = sysconf(_SC_NPROCESSORS_CONF);    // noexcept
819     auto  start           = std::chrono::steady_clock::now(); // noexcept
820
821     while (sysconf(_SC_NPROCESSORS_ONLN) < countConfigured &&
822            std::chrono::steady_clock::now() - start < std::chrono::seconds(2))
823     {
824         for (int i = 1; i < 10000; i++)
825         {
826             dummy /= i;
827         }
828     }
829
830     if (dummy < 0)
831     {
832         printf("This cannot happen, but prevents loop from being optimized away.");
833     }
834 #endif
835 }
836
837 /*! \brief Prepare the system before hardware topology detection
838  *
839  * This routine should perform any actions we want to put the system in a state
840  * where we want it to be before detecting the hardware topology. For most
841  * processors there is nothing to do, but some architectures (in particular ARM)
842  * have support for taking configured cores offline, which will make them disappear
843  * from the online processor count.
844  *
845  * This routine checks if there is a mismatch between the number of cores
846  * configured and online, and in that case we issue a small workload that
847  * attempts to wake sleeping cores before doing the actual detection.
848  *
849  * This type of mismatch can also occur for x86 or PowerPC on Linux, if SMT has only
850  * been disabled in the kernel (rather than bios). Since those cores will never
851  * come online automatically, we currently skip this test for x86 & PowerPC to
852  * avoid wasting 2 seconds. We also skip the test if there is no thread support.
853  *
854  * \note Cores will sleep relatively quickly again, so it's important to issue
855  *       the real detection code directly after this routine.
856  */
857 static void
858 hardwareTopologyPrepareDetection()
859 {
860 #if defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_CONF) && \
861     (defined(THREAD_PTHREADS) || defined(THREAD_WINDOWS))
862
863     // Modify this conditional when/if x86 or PowerPC starts to sleep some cores
864     if (!isX86 && !isPowerPC)
865     {
866         int                      countConfigured  = sysconf(_SC_NPROCESSORS_CONF);
867         std::vector<std::thread> workThreads(countConfigured);
868
869         for (auto &t : workThreads)
870         {
871             t = std::thread(spinUpCore);
872         }
873
874         for (auto &t : workThreads)
875         {
876             t.join();
877         }
878     }
879 #endif
880 }
881
882 /*! \brief Sanity check hardware topology and print some notes to log
883  *
884  *  \param mdlog            Logger.
885  *  \param hardwareTopology Reference to hardwareTopology object.
886  */
887 static void
888 hardwareTopologyDoubleCheckDetection(const gmx::MDLogger gmx_unused         &mdlog,
889                                      const gmx::HardwareTopology gmx_unused &hardwareTopology)
890 {
891 #if defined HAVE_SYSCONF && defined(_SC_NPROCESSORS_CONF)
892     if (hardwareTopology.supportLevel() < gmx::HardwareTopology::SupportLevel::LogicalProcessorCount)
893     {
894         return;
895     }
896
897     int countFromDetection = hardwareTopology.machine().logicalProcessorCount;
898     int countConfigured    = sysconf(_SC_NPROCESSORS_CONF);
899
900     /* BIOS, kernel or user actions can take physical processors
901      * offline. We already cater for the some of the cases inside the hardwareToplogy
902      * by trying to spin up cores just before we detect, but there could be other
903      * cases where it is worthwhile to hint that there might be more resources available.
904      */
905     if (countConfigured >= 0 && countConfigured != countFromDetection)
906     {
907         GMX_LOG(mdlog.info).
908             appendTextFormatted("Note: %d CPUs configured, but only %d were detected to be online.\n", countConfigured, countFromDetection);
909
910         if (isX86 && countConfigured == 2*countFromDetection)
911         {
912             GMX_LOG(mdlog.info).
913                 appendText("      X86 Hyperthreading is likely disabled; enable it for better performance.");
914         }
915         // For PowerPC (likely Power8) it is possible to set SMT to either 2,4, or 8-way hardware threads.
916         // We only warn if it is completely disabled since default performance drops with SMT8.
917         if (isPowerPC && countConfigured == 8*countFromDetection)
918         {
919             GMX_LOG(mdlog.info).
920                 appendText("      PowerPC SMT is likely disabled; enable SMT2/SMT4 for better performance.");
921         }
922     }
923 #endif
924 }
925
926
927 gmx_hw_info_t *gmx_detect_hardware(const gmx::MDLogger &mdlog, const t_commrec *cr,
928                                    gmx_bool bDetectGPUs)
929 {
930     int ret;
931
932     /* make sure no one else is doing the same thing */
933     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
934     if (ret != 0)
935     {
936         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
937     }
938
939     /* only initialize the hwinfo structure if it is not already initalized */
940     if (n_hwinfo == 0)
941     {
942         snew(hwinfo_g, 1);
943
944         hwinfo_g->cpuInfo             = new gmx::CpuInfo(gmx::CpuInfo::detect());
945
946         hardwareTopologyPrepareDetection();
947         hwinfo_g->hardwareTopology    = new gmx::HardwareTopology(gmx::HardwareTopology::detect());
948
949         // If we detected the topology on this system, double-check that it makes sense
950         if (hwinfo_g->hardwareTopology->isThisSystem())
951         {
952             hardwareTopologyDoubleCheckDetection(mdlog, *(hwinfo_g->hardwareTopology));
953         }
954
955         // TODO: Get rid of this altogether.
956         hwinfo_g->nthreads_hw_avail = hwinfo_g->hardwareTopology->machine().logicalProcessorCount;
957
958         /* detect GPUs */
959         hwinfo_g->gpu_info.n_dev            = 0;
960         hwinfo_g->gpu_info.n_dev_compatible = 0;
961         hwinfo_g->gpu_info.gpu_dev          = NULL;
962
963         /* Run the detection if the binary was compiled with GPU support
964          * and we requested detection.
965          */
966         hwinfo_g->gpu_info.bDetectGPUs =
967             (bGPUBinary && bDetectGPUs &&
968              getenv("GMX_DISABLE_GPU_DETECTION") == NULL);
969         if (hwinfo_g->gpu_info.bDetectGPUs)
970         {
971             gmx_detect_gpus(mdlog, cr);
972         }
973
974         gmx_collect_hardware_mpi(*hwinfo_g->cpuInfo);
975     }
976     /* increase the reference counter */
977     n_hwinfo++;
978
979     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
980     if (ret != 0)
981     {
982         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
983     }
984
985     return hwinfo_g;
986 }
987
988 static std::string detected_hardware_string(const gmx_hw_info_t *hwinfo,
989                                             bool                 bFullCpuInfo)
990 {
991     std::string                  s;
992
993     const gmx::CpuInfo          &cpuInfo = *hwinfo_g->cpuInfo;
994     const gmx::HardwareTopology &hwTop   = *hwinfo->hardwareTopology;
995
996     s  = gmx::formatString("\n");
997     s += gmx::formatString("Running on %d node%s with total",
998                            hwinfo->nphysicalnode,
999                            hwinfo->nphysicalnode == 1 ? "" : "s");
1000     if (hwinfo->ncore_tot > 0)
1001     {
1002         s += gmx::formatString(" %d cores,", hwinfo->ncore_tot);
1003     }
1004     s += gmx::formatString(" %d logical cores", hwinfo->nhwthread_tot);
1005     if (hwinfo->gpu_info.bDetectGPUs)
1006     {
1007         s += gmx::formatString(", %d compatible GPU%s",
1008                                hwinfo->ngpu_compatible_tot,
1009                                hwinfo->ngpu_compatible_tot == 1 ? "" : "s");
1010     }
1011     else if (bGPUBinary)
1012     {
1013         s += gmx::formatString(" (GPU detection deactivated)");
1014     }
1015     s += gmx::formatString("\n");
1016
1017     if (hwinfo->nphysicalnode > 1)
1018     {
1019         /* Print per node hardware feature counts */
1020         if (hwinfo->ncore_max > 0)
1021         {
1022             s += gmx::formatString("  Cores per node:           %2d", hwinfo->ncore_min);
1023             if (hwinfo->ncore_max > hwinfo->ncore_min)
1024             {
1025                 s += gmx::formatString(" - %2d", hwinfo->ncore_max);
1026             }
1027             s += gmx::formatString("\n");
1028         }
1029         s += gmx::formatString("  Logical cores per node:   %2d", hwinfo->nhwthread_min);
1030         if (hwinfo->nhwthread_max > hwinfo->nhwthread_min)
1031         {
1032             s += gmx::formatString(" - %2d", hwinfo->nhwthread_max);
1033         }
1034         s += gmx::formatString("\n");
1035         if (bGPUBinary)
1036         {
1037             s += gmx::formatString("  Compatible GPUs per node: %2d",
1038                                    hwinfo->ngpu_compatible_min);
1039             if (hwinfo->ngpu_compatible_max > hwinfo->ngpu_compatible_min)
1040             {
1041                 s += gmx::formatString(" - %2d", hwinfo->ngpu_compatible_max);
1042             }
1043             s += gmx::formatString("\n");
1044             if (hwinfo->ngpu_compatible_tot > 0)
1045             {
1046                 if (hwinfo->bIdenticalGPUs)
1047                 {
1048                     s += gmx::formatString("  All nodes have identical type(s) of GPUs\n");
1049                 }
1050                 else
1051                 {
1052                     /* This message will also appear with identical GPU types
1053                      * when at least one node has no GPU.
1054                      */
1055                     s += gmx::formatString("  Different nodes have different type(s) and/or order of GPUs\n");
1056                 }
1057             }
1058         }
1059     }
1060
1061 #if GMX_LIB_MPI
1062     char host[HOSTNAMELEN];
1063     int  rank;
1064
1065     gmx_gethostname(host, HOSTNAMELEN);
1066     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
1067
1068     s += gmx::formatString("Hardware detected on host %s (the node of MPI rank %d):\n",
1069                            host, rank);
1070 #else
1071     s += gmx::formatString("Hardware detected:\n");
1072 #endif
1073     s += gmx::formatString("  CPU info:\n");
1074
1075     s += gmx::formatString("    Vendor: %s\n", cpuInfo.vendorString().c_str());
1076
1077     s += gmx::formatString("    Brand:  %s\n", cpuInfo.brandString().c_str());
1078
1079     if (bFullCpuInfo)
1080     {
1081         s += gmx::formatString("    Family: %d   Model: %d   Stepping: %d\n",
1082                                cpuInfo.family(), cpuInfo.model(), cpuInfo.stepping());
1083
1084         s += gmx::formatString("    Features:");
1085         for (auto &f : cpuInfo.featureSet())
1086         {
1087             s += gmx::formatString(" %s", cpuInfo.featureString(f).c_str());;
1088         }
1089         s += gmx::formatString("\n");
1090     }
1091
1092     s += gmx::formatString("    SIMD instructions most likely to fit this hardware: %s",
1093                            gmx::simdString(static_cast<gmx::SimdType>(hwinfo->simd_suggest_min)).c_str());
1094
1095     if (hwinfo->simd_suggest_max > hwinfo->simd_suggest_min)
1096     {
1097         s += gmx::formatString(" - %s", gmx::simdString(static_cast<gmx::SimdType>(hwinfo->simd_suggest_max)).c_str());
1098     }
1099     s += gmx::formatString("\n");
1100
1101     s += gmx::formatString("    SIMD instructions selected at GROMACS compile time: %s\n",
1102                            gmx::simdString(gmx::simdCompiled()).c_str());
1103
1104     s += gmx::formatString("\n");
1105
1106     s += gmx::formatString("  Hardware topology: ");
1107     switch (hwTop.supportLevel())
1108     {
1109         case gmx::HardwareTopology::SupportLevel::None:
1110             s += gmx::formatString("None\n");
1111             break;
1112         case gmx::HardwareTopology::SupportLevel::LogicalProcessorCount:
1113             s += gmx::formatString("Only logical processor count\n");
1114             break;
1115         case gmx::HardwareTopology::SupportLevel::Basic:
1116             s += gmx::formatString("Basic\n");
1117             break;
1118         case gmx::HardwareTopology::SupportLevel::Full:
1119             s += gmx::formatString("Full\n");
1120             break;
1121         case gmx::HardwareTopology::SupportLevel::FullWithDevices:
1122             s += gmx::formatString("Full, with devices\n");
1123             break;
1124     }
1125
1126     if (!hwTop.isThisSystem())
1127     {
1128         s += gmx::formatString("  NOTE: Hardware topology cached or synthetic, not detected.\n");
1129         if (char *p = getenv("HWLOC_XMLFILE"))
1130         {
1131             s += gmx::formatString("        HWLOC_XMLFILE=%s\n", p);
1132         }
1133     }
1134
1135     if (bFullCpuInfo)
1136     {
1137         if (hwTop.supportLevel() >= gmx::HardwareTopology::SupportLevel::Basic)
1138         {
1139             s += gmx::formatString("    Sockets, cores, and logical processors:\n");
1140
1141             for (auto &socket : hwTop.machine().sockets)
1142             {
1143                 s += gmx::formatString("      Socket %2d:", socket.id);
1144                 for (auto &c : socket.cores)
1145                 {
1146                     s += gmx::formatString(" [");
1147                     for (auto &t : c.hwThreads)
1148                     {
1149                         s += gmx::formatString(" %3d", t.logicalProcessorId);
1150                     }
1151                     s += gmx::formatString("]");
1152                 }
1153                 s += gmx::formatString("\n");
1154             }
1155         }
1156         if (hwTop.supportLevel() >= gmx::HardwareTopology::SupportLevel::Full)
1157         {
1158             s += gmx::formatString("    Numa nodes:\n");
1159             for (auto &n : hwTop.machine().numa.nodes)
1160             {
1161                 s += gmx::formatString("      Node %2d (%" GMX_PRIu64 " bytes mem):", n.id, n.memory);
1162                 for (auto &l : n.logicalProcessorId)
1163                 {
1164                     s += gmx::formatString(" %3d", l);
1165                 }
1166                 s += gmx::formatString("\n");
1167             }
1168             s += gmx::formatString("      Latency:\n          ");
1169             for (std::size_t j = 0; j < hwTop.machine().numa.nodes.size(); j++)
1170             {
1171                 s += gmx::formatString(" %5d", j);
1172             }
1173             s += gmx::formatString("\n");
1174             for (std::size_t i = 0; i < hwTop.machine().numa.nodes.size(); i++)
1175             {
1176                 s += gmx::formatString("     %5d", i);
1177                 for (std::size_t j = 0; j < hwTop.machine().numa.nodes.size(); j++)
1178                 {
1179                     s += gmx::formatString(" %5.2f", hwTop.machine().numa.relativeLatency[i][j]);
1180                 }
1181                 s += gmx::formatString("\n");
1182             }
1183
1184
1185             s += gmx::formatString("    Caches:\n");
1186             for (auto &c : hwTop.machine().caches)
1187             {
1188                 s += gmx::formatString("      L%d: %" GMX_PRIu64 " bytes, linesize %d bytes, assoc. %d, shared %d ways\n",
1189                                        c.level, c.size, c.linesize, c.associativity, c.shared);
1190             }
1191         }
1192         if (hwTop.supportLevel() >= gmx::HardwareTopology::SupportLevel::FullWithDevices)
1193         {
1194             s += gmx::formatString("    PCI devices:\n");
1195             for (auto &d : hwTop.machine().devices)
1196             {
1197                 s += gmx::formatString("      %04x:%02x:%02x.%1x  Id: %04x:%04x  Class: 0x%04x  Numa: %d\n",
1198                                        d.domain, d.bus, d.dev, d.func, d.vendorId, d.deviceId, d.classId, d.numaNodeId);
1199             }
1200         }
1201     }
1202
1203     if (bGPUBinary && (hwinfo->ngpu_compatible_tot > 0 ||
1204                        hwinfo->gpu_info.n_dev > 0))
1205     {
1206         s += gmx::formatString("  GPU info:\n");
1207         s += gmx::formatString("    Number of GPUs detected: %d\n",
1208                                hwinfo->gpu_info.n_dev);
1209         if (hwinfo->gpu_info.n_dev > 0)
1210         {
1211             s += sprint_gpus(&hwinfo->gpu_info) + "\n";
1212         }
1213     }
1214     return s;
1215 }
1216
1217 void gmx_print_detected_hardware(FILE *fplog, const t_commrec *cr,
1218                                  const gmx::MDLogger &mdlog,
1219                                  const gmx_hw_info_t *hwinfo)
1220 {
1221     const gmx::CpuInfo &cpuInfo = *hwinfo_g->cpuInfo;
1222
1223     if (fplog != NULL)
1224     {
1225         std::string detected;
1226
1227         detected = detected_hardware_string(hwinfo, TRUE);
1228
1229         fprintf(fplog, "%s\n", detected.c_str());
1230     }
1231
1232     if (MULTIMASTER(cr))
1233     {
1234         std::string detected;
1235
1236         detected = detected_hardware_string(hwinfo, FALSE);
1237
1238         fprintf(stderr, "%s\n", detected.c_str());
1239     }
1240
1241     /* Check the compiled SIMD instruction set against that of the node
1242      * with the lowest SIMD level support (skip if SIMD detection did not work)
1243      */
1244     if (cpuInfo.supportLevel() >= gmx::CpuInfo::SupportLevel::Features)
1245     {
1246         gmx::simdCheck(static_cast<gmx::SimdType>(hwinfo->simd_suggest_min), fplog, MULTIMASTER(cr));
1247     }
1248
1249     /* For RDTSCP we only check on our local node and skip the MPI reduction */
1250     check_use_of_rdtscp_on_this_cpu(mdlog, cpuInfo);
1251 }
1252
1253 //! \brief Return if any GPU ID (e.g in a user-supplied string) is repeated
1254 static gmx_bool anyGpuIdIsRepeated(const gmx_gpu_opt_t *gpu_opt)
1255 {
1256     /* Loop over IDs in the string */
1257     for (int i = 0; i < gpu_opt->n_dev_use - 1; ++i)
1258     {
1259         /* Look for the ID in location i in the following part of the
1260            string */
1261         for (int j = i + 1; j < gpu_opt->n_dev_use; ++j)
1262         {
1263             if (gpu_opt->dev_use[i] == gpu_opt->dev_use[j])
1264             {
1265                 /* Same ID found in locations i and j */
1266                 return TRUE;
1267             }
1268         }
1269     }
1270
1271     return FALSE;
1272 }
1273
1274 void gmx_parse_gpu_ids(gmx_gpu_opt_t *gpu_opt)
1275 {
1276     char *env;
1277
1278     if (gpu_opt->gpu_id != NULL && !bGPUBinary)
1279     {
1280         gmx_fatal(FARGS, "GPU ID string set, but %s was compiled without GPU support!",
1281                   gmx::getProgramContext().displayName());
1282     }
1283
1284     env = getenv("GMX_GPU_ID");
1285     if (env != NULL && gpu_opt->gpu_id != NULL)
1286     {
1287         gmx_fatal(FARGS, "GMX_GPU_ID and -gpu_id can not be used at the same time");
1288     }
1289     if (env == NULL)
1290     {
1291         env = gpu_opt->gpu_id;
1292     }
1293
1294     /* parse GPU IDs if the user passed any */
1295     if (env != NULL)
1296     {
1297         /* Parse a "plain" or comma-separated GPU ID string which contains a
1298          * sequence of digits corresponding to GPU IDs; the order will
1299          * indicate the process/tMPI thread - GPU assignment. */
1300         parse_digits_from_string(env, &gpu_opt->n_dev_use, &gpu_opt->dev_use);
1301
1302         if (!gmx_multiple_gpu_per_node_supported() && 1 < gpu_opt->n_dev_use)
1303         {
1304             gmx_fatal(FARGS, "The %s implementation only supports using exactly one PP rank per node", getGpuImplementationString());
1305         }
1306         if (!gmx_gpu_sharing_supported() && anyGpuIdIsRepeated(gpu_opt))
1307         {
1308             gmx_fatal(FARGS, "The %s implementation only supports using exactly one PP rank per GPU", getGpuImplementationString());
1309         }
1310         if (gpu_opt->n_dev_use == 0)
1311         {
1312             gmx_fatal(FARGS, "Empty GPU ID string encountered.\n%s\n",
1313                       invalid_gpuid_hint);
1314         }
1315
1316         gpu_opt->bUserSet = TRUE;
1317     }
1318 }
1319
1320 void gmx_hardware_info_free(gmx_hw_info_t *hwinfo)
1321 {
1322     int ret;
1323
1324     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
1325     if (ret != 0)
1326     {
1327         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
1328     }
1329
1330     /* decrease the reference counter */
1331     n_hwinfo--;
1332
1333
1334     if (hwinfo != hwinfo_g)
1335     {
1336         gmx_incons("hwinfo < hwinfo_g");
1337     }
1338
1339     if (n_hwinfo < 0)
1340     {
1341         gmx_incons("n_hwinfo < 0");
1342     }
1343
1344     if (n_hwinfo == 0)
1345     {
1346         delete hwinfo_g->cpuInfo;
1347         delete hwinfo_g->hardwareTopology;
1348         free_gpu_info(&hwinfo_g->gpu_info);
1349         sfree(hwinfo_g);
1350     }
1351
1352     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
1353     if (ret != 0)
1354     {
1355         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
1356     }
1357 }